LCOV - code coverage report
Current view: top level - colvar - Gyration.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 147 189 77.8 %
Date: 2021-11-18 15:22:58 Functions: 10 11 90.9 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2012-2020 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 "Colvar.h"
      23             : #include "ActionRegister.h"
      24             : #include "core/PlumedMain.h"
      25             : 
      26             : #include <string>
      27             : #include <cmath>
      28             : 
      29             : using namespace std;
      30             : 
      31             : namespace PLMD {
      32             : namespace colvar {
      33             : 
      34             : //+PLUMEDOC COLVAR GYRATION
      35             : /*
      36             : Calculate the radius of gyration, or other properties related to it.
      37             : 
      38             : The different properties can be calculated and selected by the TYPE keyword:
      39             : the Radius of Gyration (RADIUS); the Trace of the Gyration Tensor (TRACE);
      40             : the Largest Principal Moment of the Gyration Tensor (GTPC_1); the middle Principal Moment of the Gyration Tensor (GTPC_2);
      41             : the Smallest Principal Moment of the Gyration Tensor (GTPC_3); the Asphericiry (ASPHERICITY); the Acylindricity (ACYLINDRICITY);
      42             : the Relative Shape Anisotropy (KAPPA2); the Smallest Principal Radius Of Gyration (GYRATION_3);
      43             : the Middle Principal Radius of Gyration (GYRATION_2); the Largest Principal Radius of Gyration (GYRATION_1).
      44             : A derivation of all these different variants can be found in \cite Vymetal:2011gv
      45             : 
      46             : The radius of gyration is calculated using:
      47             : 
      48             : \f[
      49             : s_{\rm Gyr}=\Big ( \frac{\sum_i^{n}
      50             :  m_i \vert {r}_i -{r}_{\rm COM} \vert ^2 }{\sum_i^{n} m_i} \Big)^{1/2}
      51             : \f]
      52             : 
      53             : with the position of the center of mass \f${r}_{\rm COM}\f$ given by:
      54             : 
      55             : \f[
      56             : {r}_{\rm COM}=\frac{\sum_i^{n} {r}_i\ m_i }{\sum_i^{n} m_i}
      57             : \f]
      58             : 
      59             : The radius of gyration usually makes sense when atoms used for the calculation
      60             : are all part of the same molecule.
      61             : When running with periodic boundary conditions, the atoms should be
      62             : in the proper periodic image. This is done automatically since PLUMED 2.2,
      63             : by considering the ordered list of atoms and rebuilding the broken entities using a procedure
      64             : that is equivalent to that done in \ref WHOLEMOLECULES . Notice that
      65             : rebuilding is local to this action. This is different from \ref WHOLEMOLECULES
      66             : which actually modifies the coordinates stored in PLUMED.
      67             : 
      68             : In case you want to recover the old behavior you should use the NOPBC flag.
      69             : In that case you need to take care that atoms are in the correct
      70             : periodic image.
      71             : 
      72             : 
      73             : \par Examples
      74             : 
      75             : The following input tells plumed to print the radius of gyration of the
      76             : chain containing atoms 10 to 20.
      77             : \plumedfile
      78             : GYRATION TYPE=RADIUS ATOMS=10-20 LABEL=rg
      79             : PRINT ARG=rg STRIDE=1 FILE=colvar
      80             : \endplumedfile
      81             : 
      82             : */
      83             : //+ENDPLUMEDOC
      84             : 
      85          48 : class Gyration : public Colvar {
      86             : private:
      87             :   enum CV_TYPE {RADIUS, TRACE, GTPC_1, GTPC_2, GTPC_3, ASPHERICITY, ACYLINDRICITY, KAPPA2, GYRATION_3, GYRATION_2, GYRATION_1, TOT};
      88             :   int rg_type;
      89             :   bool use_masses;
      90             :   bool nopbc;
      91             : public:
      92             :   static void registerKeywords(Keywords& keys);
      93             :   explicit Gyration(const ActionOptions&);
      94             :   virtual void calculate();
      95             : };
      96             : 
      97        7406 : PLUMED_REGISTER_ACTION(Gyration,"GYRATION")
      98             : 
      99          27 : void Gyration::registerKeywords(Keywords& keys) {
     100          27 :   Colvar::registerKeywords(keys);
     101         108 :   keys.add("atoms","ATOMS","the group of atoms that you are calculating the Gyration Tensor for");
     102         135 :   keys.add("compulsory","TYPE","RADIUS","The type of calculation relative to the Gyration Tensor you want to perform");
     103          81 :   keys.addFlag("MASS_WEIGHTED",false,"set the masses of all the atoms equal to one");
     104          27 : }
     105             : 
     106          26 : Gyration::Gyration(const ActionOptions&ao):
     107             :   PLUMED_COLVAR_INIT(ao),
     108             :   use_masses(false),
     109          28 :   nopbc(false)
     110             : {
     111             :   std::vector<AtomNumber> atoms;
     112          52 :   parseAtomList("ATOMS",atoms);
     113          25 :   if(atoms.size()==0) error("no atoms specified");
     114          50 :   parseFlag("MASS_WEIGHTED",use_masses);
     115             :   std::string Type;
     116          50 :   parse("TYPE",Type);
     117          50 :   parseFlag("NOPBC",nopbc);
     118          25 :   checkRead();
     119             : 
     120          25 :   if(Type=="RADIUS") rg_type=RADIUS;
     121          21 :   else if(Type=="TRACE") rg_type=TRACE;
     122          19 :   else if(Type=="GTPC_1") rg_type=GTPC_1;
     123          17 :   else if(Type=="GTPC_2") rg_type=GTPC_2;
     124          15 :   else if(Type=="GTPC_3") rg_type=GTPC_3;
     125          13 :   else if(Type=="ASPHERICITY") rg_type=ASPHERICITY;
     126          11 :   else if(Type=="ACYLINDRICITY") rg_type=ACYLINDRICITY;
     127           9 :   else if(Type=="KAPPA2") rg_type=KAPPA2;
     128           7 :   else if(Type=="RGYR_3") rg_type=GYRATION_3;
     129           5 :   else if(Type=="RGYR_2") rg_type=GYRATION_2;
     130           3 :   else if(Type=="RGYR_1") rg_type=GYRATION_1;
     131           2 :   else error("Unknown GYRATION type");
     132             : 
     133          24 :   switch(rg_type)
     134             :   {
     135           4 :   case RADIUS:   log.printf("  GYRATION RADIUS (Rg);"); break;
     136           2 :   case TRACE:  log.printf("  TRACE OF THE GYRATION TENSOR;"); break;
     137           2 :   case GTPC_1: log.printf("  THE LARGEST PRINCIPAL MOMENT OF THE GYRATION TENSOR (S'_1);"); break;
     138           2 :   case GTPC_2: log.printf("  THE MIDDLE PRINCIPAL MOMENT OF THE GYRATION TENSOR (S'_2);");  break;
     139           2 :   case GTPC_3: log.printf("  THE SMALLEST PRINCIPAL MOMENT OF THE GYRATION TENSOR (S'_3);"); break;
     140           2 :   case ASPHERICITY: log.printf("  THE ASPHERICITY (b');"); break;
     141           2 :   case ACYLINDRICITY: log.printf("  THE ACYLINDRICITY (c');"); break;
     142           2 :   case KAPPA2: log.printf("  THE RELATIVE SHAPE ANISOTROPY (kappa^2);"); break;
     143           2 :   case GYRATION_3: log.printf("  THE SMALLEST PRINCIPAL RADIUS OF GYRATION (r_g3);"); break;
     144           2 :   case GYRATION_2: log.printf("  THE MIDDLE PRINCIPAL RADIUS OF GYRATION (r_g2);"); break;
     145           2 :   case GYRATION_1: log.printf("  THE LARGEST PRINCIPAL RADIUS OF GYRATION (r_g1);"); break;
     146             :   }
     147          60 :   if(rg_type>TRACE) log<<"  Bibliography "<<plumed.cite("Jirí Vymetal and Jirí Vondrasek, J. Phys. Chem. A 115, 11455 (2011)");
     148          24 :   log<<"\n";
     149             : 
     150          24 :   log.printf("  atoms involved : ");
     151         444 :   for(unsigned i=0; i<atoms.size(); ++i) {
     152         132 :     if(i%25==0) log<<"\n";
     153         264 :     log.printf("%d ",atoms[i].serial());
     154             :   }
     155          24 :   log.printf("\n");
     156             : 
     157          24 :   if(nopbc) {
     158           4 :     log<<"  PBC will be ignored\n";
     159             :   } else {
     160          20 :     log<<"  broken molecules will be rebuilt assuming atoms are in the proper order\n";
     161             :   }
     162             : 
     163          24 :   addValueWithDerivatives(); setNotPeriodic();
     164          24 :   requestAtoms(atoms);
     165          24 : }
     166             : 
     167        1188 : void Gyration::calculate() {
     168             : 
     169        1188 :   if(!nopbc) makeWhole();
     170             : 
     171        1188 :   Vector com;
     172             :   double totmass = 0.;
     173        1188 :   if( use_masses ) {
     174           0 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     175           0 :       totmass+=getMass(i);
     176           0 :       com+=getMass(i)*getPosition(i);
     177             :     }
     178             :   } else {
     179        1188 :     totmass = static_cast<double>(getNumberOfAtoms());
     180       19404 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     181       18216 :       com+=getPosition(i);
     182             :     }
     183             :   }
     184        1188 :   com /= totmass;
     185             : 
     186        1188 :   double rgyr=0.;
     187        1188 :   vector<Vector> derivatives( getNumberOfAtoms() );
     188        1188 :   Tensor virial;
     189             : 
     190        1188 :   if(rg_type==RADIUS||rg_type==TRACE) {
     191         788 :     if( use_masses ) {
     192           0 :       for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     193           0 :         const Vector diff = delta( com, getPosition(i) );
     194           0 :         rgyr          += getMass(i)*diff.modulo2();
     195           0 :         derivatives[i] = diff*getMass(i);
     196           0 :         virial        -= Tensor(getPosition(i),derivatives[i]);
     197             :       }
     198             :     } else {
     199       15004 :       for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     200       14216 :         const Vector diff = delta( com, getPosition(i) );
     201        7108 :         rgyr          += diff.modulo2();
     202       14216 :         derivatives[i] = diff;
     203        7108 :         virial        -= Tensor(getPosition(i),derivatives[i]);
     204             :       }
     205             :     }
     206             :     double fact;
     207         788 :     if(rg_type==RADIUS) {
     208         658 :       rgyr = sqrt(rgyr/totmass);
     209         658 :       fact = 1./(rgyr*totmass);
     210             :     } else {
     211         130 :       rgyr = 2.*rgyr;
     212             :       fact = 4;
     213             :     }
     214         788 :     setValue(rgyr);
     215       22112 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) setAtomsDerivatives(i,fact*derivatives[i]);
     216        1576 :     setBoxDerivatives(fact*virial);
     217             :     return;
     218             :   }
     219             : 
     220             : 
     221         400 :   Tensor3d gyr_tens;
     222             :   //calculate gyration tensor
     223         400 :   if( use_masses ) {
     224           0 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     225           0 :       const Vector diff=delta( com, getPosition(i) );
     226           0 :       gyr_tens[0][0]+=getMass(i)*diff[0]*diff[0];
     227           0 :       gyr_tens[1][1]+=getMass(i)*diff[1]*diff[1];
     228           0 :       gyr_tens[2][2]+=getMass(i)*diff[2]*diff[2];
     229           0 :       gyr_tens[0][1]+=getMass(i)*diff[0]*diff[1];
     230           0 :       gyr_tens[0][2]+=getMass(i)*diff[0]*diff[2];
     231           0 :       gyr_tens[1][2]+=getMass(i)*diff[1]*diff[2];
     232             :     }
     233             :   } else {
     234        4400 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     235        4000 :       const Vector diff=delta( com, getPosition(i) );
     236        2000 :       gyr_tens[0][0]+=diff[0]*diff[0];
     237        2000 :       gyr_tens[1][1]+=diff[1]*diff[1];
     238        2000 :       gyr_tens[2][2]+=diff[2]*diff[2];
     239        2000 :       gyr_tens[0][1]+=diff[0]*diff[1];
     240        2000 :       gyr_tens[0][2]+=diff[0]*diff[2];
     241        2000 :       gyr_tens[1][2]+=diff[1]*diff[2];
     242             :     }
     243             :   }
     244             : 
     245             :   // first make the matrix symmetric
     246         400 :   gyr_tens[1][0] = gyr_tens[0][1];
     247         400 :   gyr_tens[2][0] = gyr_tens[0][2];
     248         400 :   gyr_tens[2][1] = gyr_tens[1][2];
     249         400 :   Tensor3d ttransf,transf;
     250         400 :   Vector princ_comp,prefactor;
     251             :   //diagonalize gyration tensor
     252         400 :   diagMatSym(gyr_tens, princ_comp, ttransf);
     253         400 :   transf=transpose(ttransf);
     254             :   //sort eigenvalues and eigenvectors
     255         400 :   if (princ_comp[0]<princ_comp[1]) {
     256         400 :     double tmp=princ_comp[0]; princ_comp[0]=princ_comp[1]; princ_comp[1]=tmp;
     257        1600 :     for (unsigned i=0; i<3; i++) {tmp=transf[i][0]; transf[i][0]=transf[i][1]; transf[i][1]=tmp;}
     258             :   }
     259         400 :   if (princ_comp[1]<princ_comp[2]) {
     260         400 :     double tmp=princ_comp[1]; princ_comp[1]=princ_comp[2]; princ_comp[2]=tmp;
     261        1600 :     for (unsigned i=0; i<3; i++) {tmp=transf[i][1]; transf[i][1]=transf[i][2]; transf[i][2]=tmp;}
     262             :   }
     263         400 :   if (princ_comp[0]<princ_comp[1]) {
     264         400 :     double tmp=princ_comp[0]; princ_comp[0]=princ_comp[1]; princ_comp[1]=tmp;
     265        1600 :     for (unsigned i=0; i<3; i++) {tmp=transf[i][0]; transf[i][0]=transf[i][1]; transf[i][1]=tmp;}
     266             :   }
     267             :   //calculate determinant of transformation matrix
     268             :   double det = determinant(transf);
     269             :   // trasformation matrix for rotation must have positive determinant, otherwise multiply one column by (-1)
     270         400 :   if(det<0) {
     271        1600 :     for(unsigned j=0; j<3; j++) transf[j][2]=-transf[j][2];
     272         400 :     det = -det;
     273             :   }
     274         400 :   if(fabs(det-1.)>0.0001) error("Plumed Error: Cannot diagonalize gyration tensor\n");
     275         400 :   switch(rg_type) {
     276         135 :   case GTPC_1:
     277             :   case GTPC_2:
     278             :   case GTPC_3:
     279             :   {
     280         135 :     int pc_index = rg_type-2; //index of principal component
     281         135 :     rgyr=sqrt(princ_comp[pc_index]/totmass);
     282         135 :     double rm = rgyr*totmass;
     283         135 :     if(rm>1e-6) prefactor[pc_index]=1.0/rm; //some parts of derivate
     284             :     break;
     285             :   }
     286           0 :   case GYRATION_3:        //the smallest principal radius of gyration
     287             :   {
     288           0 :     rgyr=sqrt((princ_comp[1]+princ_comp[2])/totmass);
     289           0 :     double rm = rgyr*totmass;
     290           0 :     if (rm>1e-6) {
     291           0 :       prefactor[1]=1.0/rm;
     292           0 :       prefactor[2]=1.0/rm;
     293             :     }
     294             :     break;
     295             :   }
     296         130 :   case GYRATION_2:       //the midle principal radius of gyration
     297             :   {
     298         130 :     rgyr=sqrt((princ_comp[0]+princ_comp[2])/totmass);
     299         130 :     double rm = rgyr*totmass;
     300         130 :     if (rm>1e-6) {
     301         130 :       prefactor[0]=1.0/rm;
     302         130 :       prefactor[2]=1.0/rm;
     303             :     }
     304             :     break;
     305             :   }
     306           0 :   case GYRATION_1:      //the largest principal radius of gyration
     307             :   {
     308           0 :     rgyr=sqrt((princ_comp[0]+princ_comp[1])/totmass);
     309           0 :     double rm = rgyr*totmass;
     310           0 :     if (rm>1e-6) {
     311           0 :       prefactor[0]=1.0/rm;
     312           0 :       prefactor[1]=1.0/rm;
     313             :     }
     314             :     break;
     315             :   }
     316           5 :   case ASPHERICITY:
     317             :   {
     318           5 :     rgyr=sqrt((princ_comp[0]-0.5*(princ_comp[1]+princ_comp[2]))/totmass);
     319           5 :     double rm = rgyr*totmass;
     320           5 :     if (rm>1e-6) {
     321           5 :       prefactor[0]= 1.0/rm;
     322           5 :       prefactor[1]=-0.5/rm;
     323           5 :       prefactor[2]=-0.5/rm;
     324             :     }
     325             :     break;
     326             :   }
     327           0 :   case ACYLINDRICITY:
     328             :   {
     329           0 :     rgyr=sqrt((princ_comp[1]-princ_comp[2])/totmass);
     330           0 :     double rm = rgyr*totmass;
     331           0 :     if (rm>1e-6) {  //avoid division by zero
     332           0 :       prefactor[1]= 1.0/rm;
     333           0 :       prefactor[2]=-1.0/rm;
     334             :     }
     335             :     break;
     336             :   }
     337         130 :   case KAPPA2: // relative shape anisotropy
     338             :   {
     339         130 :     double trace = princ_comp[0]+princ_comp[1]+princ_comp[2];
     340         130 :     double tmp=princ_comp[0]*princ_comp[1]+ princ_comp[1]*princ_comp[2]+ princ_comp[0]*princ_comp[2];
     341         130 :     rgyr=1.0-3*(tmp/(trace*trace));
     342         130 :     if (rgyr>1e-6) {
     343         130 :       prefactor[0]= -3*((princ_comp[1]+princ_comp[2])-2*tmp/trace)/(trace*trace) *2;
     344         130 :       prefactor[1]= -3*((princ_comp[0]+princ_comp[2])-2*tmp/trace)/(trace*trace) *2;
     345         130 :       prefactor[2]= -3*((princ_comp[0]+princ_comp[1])-2*tmp/trace)/(trace*trace) *2;
     346             :     }
     347             :     break;
     348             :   }
     349             :   }
     350             : 
     351         400 :   if(use_masses) {
     352           0 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     353           0 :       Vector tX;
     354           0 :       const Vector diff=delta( com,getPosition(i) );
     355             :       //project atomic postional vectors to diagonalized frame
     356           0 :       for(unsigned j=0; j<3; j++) tX[j]=transf[0][j]*diff[0]+transf[1][j]*diff[1]+transf[2][j]*diff[2];
     357           0 :       for(unsigned j=0; j<3; j++) derivatives[i][j]=getMass(i)*(prefactor[0]*transf[j][0]*tX[0]+
     358           0 :             prefactor[1]*transf[j][1]*tX[1]+
     359           0 :             prefactor[2]*transf[j][2]*tX[2]);
     360           0 :       setAtomsDerivatives(i,derivatives[i]);
     361             :     }
     362             :   } else {
     363        4400 :     for(unsigned i=0; i<getNumberOfAtoms(); i++) {
     364        2000 :       Vector tX;
     365        4000 :       const Vector diff=delta( com,getPosition(i) );
     366             :       //project atomic postional vectors to diagonalized frame
     367        8000 :       for(unsigned j=0; j<3; j++) tX[j]=transf[0][j]*diff[0]+transf[1][j]*diff[1]+transf[2][j]*diff[2];
     368       38000 :       for(unsigned j=0; j<3; j++) derivatives[i][j]=prefactor[0]*transf[j][0]*tX[0]+
     369       18000 :             prefactor[1]*transf[j][1]*tX[1]+
     370       12000 :             prefactor[2]*transf[j][2]*tX[2];
     371        2000 :       setAtomsDerivatives(i,derivatives[i]);
     372             :     }
     373             :   }
     374             : 
     375         400 :   setValue(rgyr);
     376         400 :   setBoxDerivativesNoPbc();
     377             : }
     378             : 
     379             : }
     380        5517 : }

Generated by: LCOV version 1.14