LCOV - code coverage report
Current view: top level - isdb - Jcoupling.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 159 171 93.0 %
Date: 2024-10-18 13:59:31 Functions: 4 5 80.0 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2016-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 "MetainferenceBase.h"
      23             : #include "core/ActionRegister.h"
      24             : #include "tools/Pbc.h"
      25             : #include "tools/Torsion.h"
      26             : 
      27             : namespace PLMD {
      28             : namespace isdb {
      29             : 
      30             : //+PLUMEDOC ISDB_COLVAR JCOUPLING
      31             : /*
      32             : Calculates 3J coupling constants for a dihedral angle.
      33             : 
      34             : The J-coupling between two atoms is given by the Karplus relation:
      35             : 
      36             : \f[
      37             : ^3J(\theta)=A\cos^2(\theta+\Delta\theta)+B\cos(\theta+\Delta\theta)+C
      38             : \f]
      39             : 
      40             : where \f$A\f$, \f$B\f$ and \f$C\f$ are the Karplus parameters and \f$\Delta\theta\f$ is an additional constant
      41             : added on to the dihedral angle \f$\theta\f$. The Karplus parameters are determined empirically and are dependent
      42             : on the type of J-coupling.
      43             : 
      44             : This collective variable computes the J-couplings for a set of atoms defining a dihedral angle. You can specify
      45             : the atoms involved using the \ref MOLINFO notation. You can also specify the experimental couplings using the
      46             :  COUPLING keywords. These will be included in the output. You must choose the type of
      47             : coupling using the type keyword, you can also supply custom Karplus parameters using TYPE=CUSTOM and the A, B, C
      48             : and SHIFT keywords. You will need to make sure you are using the correct dihedral angle:
      49             : 
      50             : - Ha-N: \f$\psi\f$
      51             : - Ha-HN: \f$\phi\f$
      52             : - N-C\f$\gamma\f$: \f$\chi_1\f$
      53             : - CO-C\f$\gamma\f$: \f$\chi_1\f$
      54             : 
      55             : J-couplings can be used to calculate a Metainference score using the internal keyword DOSCORE and all the options
      56             : of \ref METAINFERENCE .
      57             : 
      58             : \par Examples
      59             : 
      60             : In the following example we calculate the Ha-N J-coupling from a set of atoms involved in
      61             : dihedral \f$\psi\f$ angles in the peptide backbone. We also add the experimental data points and compute
      62             : the correlation and other measures and finally print the results.
      63             : 
      64             : \plumedfile
      65             : #SETTINGS MOLFILE=regtest/basic/rt32/helix.pdb
      66             : MOLINFO MOLTYPE=protein STRUCTURE=peptide.pdb
      67             : WHOLEMOLECULES ENTITY0=1-111
      68             : 
      69             : JCOUPLING ...
      70             :     TYPE=HAN
      71             :     ATOMS1=@psi-2 COUPLING1=-0.49
      72             :     ATOMS2=@psi-4 COUPLING2=-0.54
      73             :     ATOMS3=@psi-5 COUPLING3=-0.53
      74             :     ATOMS4=@psi-7 COUPLING4=-0.39
      75             :     ATOMS5=@psi-8 COUPLING5=-0.39
      76             :     LABEL=jhan
      77             : ... JCOUPLING
      78             : 
      79             : jhanst: STATS ARG=(jhan\.j-.*) PARARG=(jhan\.exp-.*)
      80             : 
      81             : PRINT ARG=jhanst.*,jhan.* FILE=COLVAR STRIDE=100
      82             : \endplumedfile
      83             : 
      84             : */
      85             : //+ENDPLUMEDOC
      86             : 
      87             : class JCoupling :
      88             :   public MetainferenceBase
      89             : {
      90             : private:
      91             :   bool pbc;
      92             :   enum { HAN, HAHN, CCG, NCG, CUSTOM };
      93             :   unsigned ncoupl_;
      94             :   double ka_;
      95             :   double kb_;
      96             :   double kc_;
      97             :   double kshift_;
      98             : 
      99             : public:
     100             :   static void registerKeywords(Keywords& keys);
     101             :   explicit JCoupling(const ActionOptions&);
     102             :   void calculate() override;
     103             :   void update() override;
     104             : };
     105             : 
     106             : PLUMED_REGISTER_ACTION(JCoupling, "JCOUPLING")
     107             : 
     108           8 : void JCoupling::registerKeywords(Keywords& keys) {
     109           8 :   MetainferenceBase::registerKeywords(keys);
     110          16 :   keys.addFlag("NOPBC",false,"ignore the periodic boundary conditions when calculating distances");
     111          16 :   keys.add("numbered", "ATOMS", "the 4 atoms involved in each of the bonds for which you wish to calculate the J-coupling. "
     112             :            "Keywords like ATOMS1, ATOMS2, ATOMS3,... should be listed and one J-coupling will be "
     113             :            "calculated for each ATOMS keyword you specify.");
     114          16 :   keys.reset_style("ATOMS", "atoms");
     115          16 :   keys.add("compulsory", "TYPE", "Type of J-coupling to compute (HAN,HAHN,CCG,NCG,CUSTOM)");
     116          16 :   keys.add("optional", "A", "Karplus parameter A");
     117          16 :   keys.add("optional", "B", "Karplus parameter B");
     118          16 :   keys.add("optional", "C", "Karplus parameter C");
     119          16 :   keys.add("optional", "SHIFT", "Angle shift in radians");
     120          16 :   keys.add("numbered", "COUPLING", "Add an experimental value for each coupling");
     121          16 :   keys.addOutputComponent("j", "default","scalar", "the calculated J-coupling");
     122          16 :   keys.addOutputComponent("exp", "COUPLING","scalar", "the experimental J-coupling");
     123           8 : }
     124             : 
     125           6 : JCoupling::JCoupling(const ActionOptions&ao):
     126             :   PLUMED_METAINF_INIT(ao),
     127           6 :   pbc(true)
     128             : {
     129           6 :   bool nopbc = !pbc;
     130           6 :   parseFlag("NOPBC", nopbc);
     131           6 :   pbc =! nopbc;
     132             : 
     133             :   // Read in the atoms
     134             :   std::vector<AtomNumber> t, atoms;
     135           6 :   for (int i = 1; ; ++i) {
     136          68 :     parseAtomList("ATOMS", i, t );
     137          34 :     if (t.empty()) {
     138             :       break;
     139             :     }
     140             : 
     141          28 :     if (t.size() != 4) {
     142             :       std::string ss;
     143           0 :       Tools::convert(i, ss);
     144           0 :       error("ATOMS" + ss + " keyword has the wrong number of atoms");
     145             :     }
     146             : 
     147             :     // This makes the distance calculation easier later on (see Torsion implementation)
     148          28 :     atoms.push_back(t[0]);
     149          28 :     atoms.push_back(t[1]);
     150          28 :     atoms.push_back(t[1]);
     151          28 :     atoms.push_back(t[2]);
     152          28 :     atoms.push_back(t[2]);
     153          28 :     atoms.push_back(t[3]);
     154          28 :     t.resize(0);
     155          28 :   }
     156             : 
     157             :   // We now have 6 atoms per datapoint
     158           6 :   ncoupl_ = atoms.size()/6;
     159             : 
     160             :   // Parse J-Coupling type, this will determine the Karplus parameters
     161             :   unsigned jtype_ = CUSTOM;
     162             :   std::string string_type;
     163          12 :   parse("TYPE", string_type);
     164           6 :   if(string_type == "HAN") {
     165             :     jtype_ = HAN;
     166           5 :   } else if(string_type == "HAHN") {
     167             :     jtype_ = HAHN;
     168           2 :   } else if(string_type == "CCG") {
     169             :     jtype_ = CCG;
     170           1 :   } else if(string_type == "NCG") {
     171             :     jtype_ = NCG;
     172           0 :   } else if(string_type == "CUSTOM") {
     173             :     jtype_ = CUSTOM;
     174             :   } else {
     175           0 :     error("Unknown J-coupling type!");
     176             :   }
     177             : 
     178             :   // Optionally add an experimental value (like with RDCs)
     179             :   std::vector<double> coupl;
     180           6 :   coupl.resize( ncoupl_ );
     181             :   unsigned ntarget=0;
     182          13 :   for(unsigned i=0; i<ncoupl_; ++i) {
     183          24 :     if( !parseNumbered( "COUPLING", i+1, coupl[i] ) ) break;
     184           7 :     ntarget++;
     185             :   }
     186             :   bool addcoupling=false;
     187           6 :   if(ntarget!=ncoupl_ && ntarget!=0) error("found wrong number of COUPLING values");
     188           6 :   if(ntarget==ncoupl_) addcoupling=true;
     189           6 :   if(getDoScore()&&!addcoupling) error("with DOSCORE you need to set the COUPLING values");
     190             : 
     191             :   // For custom types we allow use of custom Karplus parameters
     192           6 :   if (jtype_ == CUSTOM) {
     193           0 :     parse("A", ka_);
     194           0 :     parse("B", kb_);
     195           0 :     parse("C", kc_);
     196           0 :     parse("SHIFT", kshift_);
     197             :   }
     198             : 
     199           6 :   log << "  Bibliography ";
     200          12 :   log<<plumed.cite("Bonomi, Camilloni, Bioinformatics, 33, 3999 (2017)");
     201             : 
     202             :   // Set Karplus parameters
     203           6 :   switch (jtype_) {
     204           1 :   case HAN:
     205           1 :     ka_ = -0.88;
     206           1 :     kb_ = -0.61;
     207           1 :     kc_ = -0.27;
     208           1 :     kshift_ = pi / 3.0;
     209           2 :     log << plumed.cite("Wang A C, Bax A, J. Am. Chem. Soc. 117, 1810 (1995)");
     210           1 :     log<<"\n";
     211           1 :     log.printf("  J-coupling type is HAN, with A: %f, B: %f, C: %f, angle shift: %f\n", ka_, kb_, kc_, kshift_);
     212             :     break;
     213           3 :   case HAHN:
     214           3 :     ka_ = 7.09;
     215           3 :     kb_ = -1.42;
     216           3 :     kc_ = 1.55;
     217           3 :     kshift_ = -pi / 3.0;
     218           6 :     log << plumed.cite("Hu J-S, Bax A, J. Am. Chem. Soc. 119, 6360 (1997)");
     219           3 :     log<<"\n";
     220           3 :     log.printf("  J-coupling type is HAHN, with A: %f, B: %f, C: %f, angle shift: %f\n", ka_, kb_, kc_, kshift_);
     221             :     break;
     222           1 :   case CCG:
     223           1 :     ka_ = 2.31;
     224           1 :     kb_ = -0.87;
     225           1 :     kc_ = 0.55;
     226           1 :     kshift_ = (2.0 * pi) / 3.0;
     227           2 :     log << plumed.cite("Perez C, Löhr F, Rüterjans H, Schmidt J, J. Am. Chem. Soc. 123, 7081 (2001)");
     228           1 :     log<<"\n";
     229           1 :     log.printf("  J-coupling type is CCG, with A: %f, B: %f, C: %f, angle shift: %f\n", ka_, kb_, kc_, kshift_);
     230             :     break;
     231           1 :   case NCG:
     232           1 :     ka_ = 1.29;
     233           1 :     kb_ = -0.49;
     234           1 :     kc_ = 0.37;
     235           1 :     kshift_ = 0.0;
     236           2 :     log << plumed.cite("Perez C, Löhr F, Rüterjans H, Schmidt J, J. Am. Chem. Soc. 123, 7081 (2001)");
     237           1 :     log<<"\n";
     238           1 :     log.printf("  J-coupling type is NCG, with A: %f, B: %f, C: %f, angle shift: %f\n", ka_, kb_, kc_, kshift_);
     239             :     break;
     240           0 :   case CUSTOM:
     241           0 :     log<<"\n";
     242           0 :     log.printf("  J-coupling type is custom, with A: %f, B: %f, C: %f, angle shift: %f\n", ka_, kb_, kc_, kshift_);
     243             :     break;
     244             :   }
     245             : 
     246          34 :   for (unsigned i = 0; i < ncoupl_; ++i) {
     247          28 :     log.printf("  The %uth J-Coupling is calculated from atoms : %d %d %d %d.",
     248          28 :                i+1, atoms[6*i].serial(), atoms[6*i+1].serial(), atoms[6*i+3].serial(), atoms[6*i+5].serial());
     249          28 :     if (addcoupling) {
     250           7 :       log.printf(" Experimental J-Coupling is %f.", coupl[i]);
     251             :     }
     252          28 :     log.printf("\n");
     253             :   }
     254             : 
     255           6 :   if (pbc) {
     256           0 :     log.printf("  using periodic boundary conditions\n");
     257             :   } else {
     258           6 :     log.printf("  without periodic boundary conditions\n");
     259             :   }
     260             : 
     261           6 :   if(!getDoScore()) {
     262          26 :     for (unsigned i = 0; i < ncoupl_; i++) {
     263          21 :       std::string num; Tools::convert(i, num);
     264          42 :       addComponentWithDerivatives("j-" + num);
     265          42 :       componentIsNotPeriodic("j-" + num);
     266             :     }
     267             :   } else {
     268           8 :     for (unsigned i = 0; i < ncoupl_; i++) {
     269           7 :       std::string num; Tools::convert(i, num);
     270          14 :       addComponent("j-" + num);
     271          14 :       componentIsNotPeriodic("j-" + num);
     272             :     }
     273             :   }
     274             : 
     275           6 :   if (addcoupling||getDoScore()) {
     276           8 :     for (unsigned i = 0; i < ncoupl_; i++) {
     277           7 :       std::string num; Tools::convert(i, num);
     278          14 :       addComponent("exp-" + num);
     279           7 :       componentIsNotPeriodic("exp-" + num);
     280           7 :       Value* comp = getPntrToComponent("exp-" + num);
     281           7 :       comp->set(coupl[i]);
     282             :     }
     283             :   }
     284             : 
     285           6 :   requestAtoms(atoms, false);
     286           6 :   if(getDoScore()) {
     287           1 :     setParameters(coupl);
     288           1 :     Initialise(ncoupl_);
     289             :   }
     290           6 :   setDerivatives();
     291           6 :   checkRead();
     292           6 : }
     293             : 
     294          16 : void JCoupling::calculate()
     295             : {
     296          16 :   if (pbc) makeWhole();
     297          16 :   std::vector<Vector> deriv(ncoupl_*6);
     298          16 :   std::vector<double> j(ncoupl_,0.);
     299             : 
     300          16 :   #pragma omp parallel num_threads(OpenMP::getNumThreads())
     301             :   {
     302             :     #pragma omp for
     303             :     // Loop through atoms, with steps of 6 atoms (one iteration per datapoint)
     304             :     for (unsigned r=0; r<ncoupl_; r++) {
     305             :       // Index is the datapoint index
     306             :       unsigned a0 = 6*r;
     307             : 
     308             :       // 6 atoms -> 3 vectors
     309             :       Vector d0 = delta(getPosition(a0+1), getPosition(a0));
     310             :       Vector d1 = delta(getPosition(a0+3), getPosition(a0+2));
     311             :       Vector d2 = delta(getPosition(a0+5), getPosition(a0+4));
     312             : 
     313             :       // Calculate dihedral with 3 vectors, get the derivatives
     314             :       Vector dd0, dd1, dd2;
     315             :       PLMD::Torsion t;
     316             :       double torsion = t.compute(d0, d1, d2, dd0, dd1, dd2);
     317             : 
     318             :       // Calculate the Karplus relation and its derivative
     319             :       double theta = torsion + kshift_;
     320             :       double cos_theta = std::cos(theta);
     321             :       double sin_theta = std::sin(theta);
     322             :       j[r] = ka_*cos_theta*cos_theta + kb_*cos_theta + kc_;
     323             :       double derj = -2.*ka_*sin_theta*cos_theta - kb_*sin_theta;
     324             : 
     325             :       dd0 *= derj;
     326             :       dd1 *= derj;
     327             :       dd2 *= derj;
     328             : 
     329             :       if(getDoScore()) setCalcData(r, j[r]);
     330             :       deriv[a0] =  dd0;
     331             :       deriv[a0+1] = -dd0;
     332             :       deriv[a0+2] =  dd1;
     333             :       deriv[a0+3] = -dd1;
     334             :       deriv[a0+4] =  dd2;
     335             :       deriv[a0+5] = -dd2;
     336             :     }
     337             :   }
     338             : 
     339          16 :   if(getDoScore()) {
     340             :     /* Metainference */
     341           6 :     double score = getScore();
     342           6 :     setScore(score);
     343             : 
     344             :     /* calculate final derivatives */
     345           6 :     Tensor virial;
     346           6 :     Value* val=getPntrToComponent("score");
     347          48 :     for (unsigned r=0; r<ncoupl_; r++) {
     348          42 :       const unsigned a0 = 6*r;
     349          42 :       setAtomsDerivatives(val, a0, deriv[a0]*getMetaDer(r));
     350          42 :       setAtomsDerivatives(val, a0+1, deriv[a0+1]*getMetaDer(r));
     351          42 :       setAtomsDerivatives(val, a0+2, deriv[a0+2]*getMetaDer(r));
     352          42 :       setAtomsDerivatives(val, a0+3, deriv[a0+3]*getMetaDer(r));
     353          42 :       setAtomsDerivatives(val, a0+4, deriv[a0+4]*getMetaDer(r));
     354          42 :       setAtomsDerivatives(val, a0+5, deriv[a0+5]*getMetaDer(r));
     355          42 :       virial-=Tensor(getPosition(a0), deriv[a0]*getMetaDer(r));
     356          42 :       virial-=Tensor(getPosition(a0+1), deriv[a0+1]*getMetaDer(r));
     357          42 :       virial-=Tensor(getPosition(a0+2), deriv[a0+2]*getMetaDer(r));
     358          42 :       virial-=Tensor(getPosition(a0+3), deriv[a0+3]*getMetaDer(r));
     359          42 :       virial-=Tensor(getPosition(a0+4), deriv[a0+4]*getMetaDer(r));
     360          42 :       virial-=Tensor(getPosition(a0+5), deriv[a0+5]*getMetaDer(r));
     361             :     }
     362           6 :     setBoxDerivatives(val, virial);
     363             :   } else {
     364          66 :     for (unsigned r=0; r<ncoupl_; r++) {
     365          56 :       const unsigned a0 = 6*r;
     366          56 :       std::string num; Tools::convert(r,num);
     367          56 :       Value* val=getPntrToComponent("j-"+num);
     368          56 :       val->set(j[r]);
     369          56 :       setAtomsDerivatives(val, a0, deriv[a0]);
     370          56 :       setAtomsDerivatives(val, a0+1, deriv[a0+1]);
     371          56 :       setAtomsDerivatives(val, a0+2, deriv[a0+2]);
     372          56 :       setAtomsDerivatives(val, a0+3, deriv[a0+3]);
     373          56 :       setAtomsDerivatives(val, a0+4, deriv[a0+4]);
     374          56 :       setAtomsDerivatives(val, a0+5, deriv[a0+5]);
     375          56 :       Tensor virial;
     376          56 :       virial-=Tensor(getPosition(a0), deriv[a0]);
     377          56 :       virial-=Tensor(getPosition(a0+1), deriv[a0+1]);
     378          56 :       virial-=Tensor(getPosition(a0+2), deriv[a0+2]);
     379          56 :       virial-=Tensor(getPosition(a0+3), deriv[a0+3]);
     380          56 :       virial-=Tensor(getPosition(a0+4), deriv[a0+4]);
     381          56 :       virial-=Tensor(getPosition(a0+5), deriv[a0+5]);
     382          56 :       setBoxDerivatives(val, virial);
     383             :     }
     384             :   }
     385          16 : }
     386             : 
     387          16 : void JCoupling::update() {
     388             :   // write status file
     389          16 :   if(getWstride()>0&& (getStep()%getWstride()==0 || getCPT()) ) writeStatus();
     390          16 : }
     391             : 
     392             : }
     393             : }

Generated by: LCOV version 1.16