LCOV - code coverage report
Current view: top level - isdb - RDC.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 180 197 91.4 %
Date: 2025-03-25 09:33:27 Functions: 5 6 83.3 %

          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 "MetainferenceBase.h"
      23             : #include "core/ActionRegister.h"
      24             : #include "tools/Pbc.h"
      25             : 
      26             : #ifdef __PLUMED_HAS_GSL
      27             : #include <gsl/gsl_vector.h>
      28             : #include <gsl/gsl_matrix.h>
      29             : #include <gsl/gsl_linalg.h>
      30             : #include <gsl/gsl_blas.h>
      31             : #endif
      32             : 
      33             : namespace PLMD {
      34             : namespace isdb {
      35             : 
      36             : //+PLUMEDOC ISDB_COLVAR RDC
      37             : /*
      38             : Calculates the (Residual) Dipolar Coupling between two atoms.
      39             : 
      40             : The Dipolar Coupling between two nuclei depends on the \f$\theta\f$ angle between
      41             : the inter-nuclear vector and the external magnetic field.
      42             : 
      43             : \f[
      44             : D=D_{max}0.5(3\cos^2(\theta)-1)
      45             : \f]
      46             : 
      47             : where
      48             : 
      49             : \f[
      50             : D_{max}=-\mu_0\gamma_1\gamma_2h/(8\pi^3r^3)
      51             : \f]
      52             : 
      53             : that is the maximal value of the dipolar coupling for the two nuclear spins with gyromagnetic ratio \f$\gamma\f$.
      54             : \f$\mu\f$ is the magnetic constant and h is the Planck constant.
      55             : 
      56             : Common Gyromagnetic Ratios (C.G.S)
      57             : - H(1) 26.7513
      58             : - C(13) 6.7261
      59             : - N(15) -2.7116
      60             : and their products (this is what is given in input using the keyword GYROM)
      61             : - N-H -72.5388
      62             : - C-H 179.9319
      63             : - C-N -18.2385
      64             : - C-C 45.2404
      65             : 
      66             : In isotropic media DCs average to zero because of the rotational
      67             : averaging, but when the rotational symmetry is broken, either through the introduction of an alignment medium or for molecules
      68             : with highly anisotropic paramagnetic susceptibility, then the average of the DCs is not zero and it is possible to measure a Residual Dipolar Coupling (RDCs).
      69             : 
      70             : This collective variable calculates the Dipolar Coupling for a set of couple of atoms using
      71             : the above definition.
      72             : 
      73             : In a standard MD simulation the average over time of the DC should then be zero. If one wants to model the meaning of a set of measured RDCs it is possible to try to solve the following problem: "what is the distribution of structures and orientations that reproduce the measured RDCs".
      74             : 
      75             : This collective variable can then be use to break the rotational symmetry of a simulation by imposing that the average of the DCs over the conformational ensemble must be equal to the measured RDCs \cite Camilloni:2015ka . Since measured RDCs are also a function of the fraction of aligned molecules in the sample it is better to compare them modulo a constant or looking at the correlation.
      76             : 
      77             : Alternatively if the molecule is rigid it is possible to use the experimental data to calculate the alignment tensor and the use that to back calculate the RDCs, this is what is usually call the Single Value Decomposition approach. In this case the code rely on the
      78             : a set of function from the GNU Scientific Library (GSL). (With SVD forces are not currently implemented).
      79             : 
      80             : Replica-Averaged simulations can be performed using RDCs, \ref ENSEMBLE, \ref STATS and \ref RESTRAINT .
      81             : \ref METAINFERENCE can be activated using DOSCORE and the other relevant keywords.
      82             : 
      83             : Additional material and examples can be also found in the tutorial \ref isdb-1
      84             : 
      85             : \par Examples
      86             : In the following example five N-H RDCs are defined and averaged over multiple replicas,
      87             : their correlation is then calculated with respect to a set of experimental data and restrained.
      88             : In addition, and only for analysis purposes, the same RDCs each single conformation are calculated
      89             : using a Single Value Decomposition algorithm, then averaged and again compared with the experimental data.
      90             : 
      91             : \plumedfile
      92             : #SETTINGS NREPLICAS=2
      93             : RDC ...
      94             : GYROM=-72.5388
      95             : SCALE=0.001
      96             : ATOMS1=20,21
      97             : ATOMS2=37,38
      98             : ATOMS3=56,57
      99             : ATOMS4=76,77
     100             : ATOMS5=92,93
     101             : LABEL=nh
     102             : ... RDC
     103             : 
     104             : erdc: ENSEMBLE ARG=nh.*
     105             : 
     106             : st: STATS ARG=erdc.* PARAMETERS=8.17,-8.271,-10.489,-9.871,-9.152
     107             : 
     108             : rdce: RESTRAINT ARG=st.corr KAPPA=0. SLOPE=-25000.0 AT=1.
     109             : 
     110             : RDC ...
     111             : GYROM=-72.5388
     112             : SVD
     113             : ATOMS1=20,21 COUPLING1=8.17
     114             : ATOMS2=37,38 COUPLING2=-8.271
     115             : ATOMS3=56,57 COUPLING3=-10.489
     116             : ATOMS4=76,77 COUPLING4=-9.871
     117             : ATOMS5=92,93 COUPLING5=-9.152
     118             : LABEL=svd
     119             : ... RDC
     120             : 
     121             : esvd: ENSEMBLE ARG=(svd\.rdc-.*)
     122             : 
     123             : st_svd: STATS ARG=esvd.* PARAMETERS=8.17,-8.271,-10.489,-9.871,-9.152
     124             : 
     125             : PRINT ARG=st.corr,st_svd.corr,rdce.bias FILE=colvar
     126             : \endplumedfile
     127             : 
     128             : */
     129             : //+ENDPLUMEDOC
     130             : 
     131             : //+PLUMEDOC ISDB_COLVAR PCS
     132             : /*
     133             : Calculates the Pseudo-contact shift of a nucleus determined by the presence of a metal ion susceptible to anisotropic magnetization.
     134             : 
     135             : The PCS of an atomic nucleus depends on the \f$\theta\f$ angle between the vector from the spin-label to the nucleus
     136             :  and the external magnetic field and the module of the vector itself \cite Camilloni:2015jf . While in principle the averaging
     137             : resulting from the tumbling should remove the pseudo-contact shift, in presence of the NMR magnetic field the magnetically anisotropic molecule bound to system will break the rotational symmetry does resulting in measurable values for the PCS and RDC.
     138             : 
     139             : PCS values can also be calculated using a Single Value Decomposition approach, in this case the code rely on the
     140             : a set of function from the GNU Scientific Library (GSL). (With SVD forces are not currently implemented).
     141             : 
     142             : Replica-Averaged simulations can be performed using PCS values, \ref ENSEMBLE, \ref STATS and \ref RESTRAINT .
     143             : Metainference simulations can be performed with this CV and \ref METAINFERENCE .
     144             : 
     145             : \par Examples
     146             : 
     147             : In the following example five PCS values are defined and their correlation with
     148             : respect to a set of experimental data is calculated and restrained. In addition,
     149             : and only for analysis purposes, the same PCS values are calculated using a Single Value
     150             : Decomposition algorithm.
     151             : 
     152             : \plumedfile
     153             : PCS ...
     154             : ATOMS1=20,21
     155             : ATOMS2=20,38
     156             : ATOMS3=20,57
     157             : ATOMS4=20,77
     158             : ATOMS5=20,93
     159             : LABEL=nh
     160             : ... PCS
     161             : 
     162             : enh: ENSEMBLE ARG=nh.*
     163             : 
     164             : st: STATS ARG=enh.* PARAMETERS=8.17,-8.271,-10.489,-9.871,-9.152
     165             : 
     166             : pcse: RESTRAINT ARG=st.corr KAPPA=0. SLOPE=-25000.0 AT=1.
     167             : 
     168             : PRINT ARG=st.corr,pcse.bias FILE=colvar
     169             : \endplumedfile
     170             : 
     171             : */
     172             : //+ENDPLUMEDOC
     173             : 
     174             : class RDC :
     175             :   public MetainferenceBase {
     176             : private:
     177             :   double         Const;
     178             :   double         mu_s;
     179             :   double         scale;
     180             :   std::vector<double> coupl;
     181             :   bool           svd;
     182             :   bool           pbc;
     183             : 
     184             : #ifdef __PLUMED_HAS_GSL
     185             : /// Auxiliary class to delete a gsl_vector.
     186             : /// If used somewhere else we can move it.
     187             :   struct gsl_vector_deleter {
     188             :     void operator()(gsl_vector* p) {
     189          25 :       gsl_vector_free(p);
     190          25 :     }
     191             :   };
     192             : 
     193             : /// unique_ptr to a gsl_vector.
     194             : /// Gets deleted when going out of scope.
     195             :   typedef std::unique_ptr<gsl_vector,gsl_vector_deleter> gsl_vector_unique_ptr;
     196             : 
     197             : /// Auxiliary class to delete a gsl_matrix.
     198             : /// If used somewhere else we can move it.
     199             :   struct gsl_matrix_deleter {
     200             :     void operator()(gsl_matrix* p) {
     201          15 :       gsl_matrix_free(p);
     202          15 :     }
     203             :   };
     204             : 
     205             : /// unique_ptr to a gsl_matrix.
     206             : /// Gets deleted when going out of scope.
     207             :   typedef std::unique_ptr<gsl_matrix,gsl_matrix_deleter> gsl_matrix_unique_ptr;
     208             : #endif
     209             : 
     210             : 
     211             :   void do_svd();
     212             : public:
     213             :   explicit RDC(const ActionOptions&);
     214             :   static void registerKeywords( Keywords& keys );
     215             :   void calculate() override;
     216             :   void update() override;
     217             : };
     218             : 
     219             : PLUMED_REGISTER_ACTION(RDC,"RDC")
     220             : PLUMED_REGISTER_ACTION(RDC,"PCS")
     221             : 
     222          33 : void RDC::registerKeywords( Keywords& keys ) {
     223          33 :   MetainferenceBase::registerKeywords(keys);
     224          33 :   keys.addFlag("NOPBC",false,"ignore the periodic boundary conditions when calculating distances");
     225          33 :   keys.add("numbered","ATOMS","the couple of atoms involved in each of the bonds for which you wish to calculate the RDC. "
     226             :            "Keywords like ATOMS1, ATOMS2, ATOMS3,... should be listed and one dipolar coupling will be "
     227             :            "calculated for each ATOMS keyword you specify.");
     228          66 :   keys.reset_style("ATOMS","atoms");
     229          33 :   keys.add("compulsory","GYROM","1.","Add the product of the gyromagnetic constants for the bond. ");
     230          33 :   keys.add("compulsory","SCALE","1.","Add the scaling factor to take into account concentration and other effects. ");
     231          33 :   keys.addFlag("SVD",false,"Set to TRUE if you want to back calculate using Single Value Decomposition (need GSL at compilation time).");
     232          33 :   keys.add("numbered","COUPLING","Add an experimental value for each coupling (needed by SVD and useful for STATS).");
     233          66 :   keys.addOutputComponent("rdc","default","scalar","the calculated # RDC");
     234          66 :   keys.addOutputComponent("exp","SVD/COUPLING","scalar","the experimental # RDC");
     235          66 :   keys.addOutputComponent("Sxx","SVD","scalar","Tensor component");
     236          66 :   keys.addOutputComponent("Syy","SVD","scalar","Tensor component");
     237          66 :   keys.addOutputComponent("Szz","SVD","scalar","Tensor component");
     238          66 :   keys.addOutputComponent("Sxy","SVD","scalar","Tensor component");
     239          66 :   keys.addOutputComponent("Sxz","SVD","scalar","Tensor component");
     240          66 :   keys.addOutputComponent("Syz","SVD","scalar","Tensor component");
     241          33 : }
     242             : 
     243          29 : RDC::RDC(const ActionOptions&ao):
     244             :   PLUMED_METAINF_INIT(ao),
     245          29 :   Const(1.),
     246          29 :   mu_s(1.),
     247          29 :   scale(1.),
     248          29 :   pbc(true) {
     249          29 :   bool nopbc=!pbc;
     250          29 :   parseFlag("NOPBC",nopbc);
     251          29 :   pbc=!nopbc;
     252             : 
     253             :   const double RDCConst = 0.3356806;
     254             :   const double PCSConst = 1.0;
     255             : 
     256          29 :   if( getName().find("RDC")!=std::string::npos) {
     257          29 :     Const *= RDCConst;
     258           0 :   } else if( getName().find("PCS")!=std::string::npos) {
     259             :     Const *= PCSConst;
     260             :   }
     261             : 
     262             :   // Read in the atoms
     263             :   std::vector<AtomNumber> t, atoms;
     264          29 :   for(int i=1;; ++i ) {
     265         292 :     parseAtomList("ATOMS", i, t );
     266         146 :     if( t.empty() ) {
     267             :       break;
     268             :     }
     269         117 :     if( t.size()!=2 ) {
     270             :       std::string ss;
     271           0 :       Tools::convert(i,ss);
     272           0 :       error("ATOMS" + ss + " keyword has the wrong number of atoms");
     273             :     }
     274         117 :     atoms.push_back(t[0]);
     275         117 :     atoms.push_back(t[1]);
     276         117 :     t.resize(0);
     277         117 :   }
     278             : 
     279          29 :   const unsigned ndata = atoms.size()/2;
     280             : 
     281             :   // Read in GYROMAGNETIC constants
     282          29 :   parse("GYROM", mu_s);
     283          29 :   if(mu_s==0.) {
     284           0 :     error("GYROM cannot be 0");
     285             :   }
     286             : 
     287             :   // Read in SCALING factors
     288          29 :   parse("SCALE", scale);
     289          29 :   if(scale==0.) {
     290           0 :     error("SCALE cannot be 0");
     291             :   }
     292             : 
     293          29 :   svd=false;
     294          29 :   parseFlag("SVD",svd);
     295             : #ifndef __PLUMED_HAS_GSL
     296             :   if(svd) {
     297             :     error("You CANNOT use SVD without GSL. Recompile PLUMED with GSL!\n");
     298             :   }
     299             : #endif
     300          29 :   if(svd&&getDoScore()) {
     301           0 :     error("It is not possible to use SVD and METAINFERENCE together");
     302             :   }
     303             : 
     304             :   // Optionally add an experimental value
     305          29 :   coupl.resize( ndata );
     306             :   unsigned ntarget=0;
     307          82 :   for(unsigned i=0; i<ndata; ++i) {
     308         138 :     if( !parseNumbered( "COUPLING", i+1, coupl[i] ) ) {
     309             :       break;
     310             :     }
     311          53 :     ntarget++;
     312             :   }
     313             :   bool addexp=false;
     314          29 :   if(ntarget!=ndata && ntarget!=0) {
     315           0 :     error("found wrong number of COUPLING values");
     316             :   }
     317          29 :   if(ntarget==ndata) {
     318             :     addexp=true;
     319             :   }
     320          29 :   if(getDoScore()&&!addexp) {
     321           0 :     error("with DOSCORE you need to set the COUPLING values");
     322             :   }
     323          29 :   if(svd&&!addexp) {
     324           0 :     error("with SVD you need to set the COUPLING values");
     325             :   }
     326             : 
     327             : 
     328             :   // Output details of all contacts
     329          29 :   log.printf("  Gyromagnetic moment is %f. Scaling factor is %f.",mu_s,scale);
     330         146 :   for(unsigned i=0; i<ndata; ++i) {
     331         117 :     log.printf("  The %uth Bond Dipolar Coupling is calculated from atoms : %d %d.", i+1, atoms[2*i].serial(), atoms[2*i+1].serial());
     332         117 :     if(addexp) {
     333          53 :       log.printf(" Experimental coupling is %f.", coupl[i]);
     334             :     }
     335         117 :     log.printf("\n");
     336             :   }
     337             : 
     338          29 :   log<<"  Bibliography "
     339          58 :      <<plumed.cite("Camilloni C, Vendruscolo M, J. Phys. Chem. B 119, 653 (2015)")
     340          87 :      <<plumed.cite("Camilloni C, Vendruscolo M, Biochemistry 54, 7470 (2015)");
     341          58 :   log<< plumed.cite("Bonomi, Camilloni, Bioinformatics, 33, 3999 (2017)") << "\n";
     342             : 
     343             : 
     344          29 :   if(!getDoScore()&&!svd) {
     345          80 :     for(unsigned i=0; i<ndata; i++) {
     346             :       std::string num;
     347          64 :       Tools::convert(i,num);
     348         128 :       addComponentWithDerivatives("rdc-"+num);
     349         128 :       componentIsNotPeriodic("rdc-"+num);
     350             :     }
     351          16 :     if(addexp) {
     352           0 :       for(unsigned i=0; i<ndata; i++) {
     353             :         std::string num;
     354           0 :         Tools::convert(i,num);
     355           0 :         addComponent("exp-"+num);
     356           0 :         componentIsNotPeriodic("exp-"+num);
     357           0 :         Value* comp=getPntrToComponent("exp-"+num);
     358           0 :         comp->set(coupl[i]);
     359             :       }
     360             :     }
     361             :   } else {
     362          66 :     for(unsigned i=0; i<ndata; i++) {
     363             :       std::string num;
     364          53 :       Tools::convert(i,num);
     365         106 :       addComponentWithDerivatives("rdc-"+num);
     366         106 :       componentIsNotPeriodic("rdc-"+num);
     367             :     }
     368          66 :     for(unsigned i=0; i<ndata; i++) {
     369             :       std::string num;
     370          53 :       Tools::convert(i,num);
     371         106 :       addComponent("exp-"+num);
     372          53 :       componentIsNotPeriodic("exp-"+num);
     373          53 :       Value* comp=getPntrToComponent("exp-"+num);
     374          53 :       comp->set(coupl[i]);
     375             :     }
     376             :   }
     377             : 
     378          29 :   if(svd) {
     379           2 :     addComponent("Sxx");
     380           2 :     componentIsNotPeriodic("Sxx");
     381           2 :     addComponent("Syy");
     382           2 :     componentIsNotPeriodic("Syy");
     383           2 :     addComponent("Szz");
     384           2 :     componentIsNotPeriodic("Szz");
     385           2 :     addComponent("Sxy");
     386           2 :     componentIsNotPeriodic("Sxy");
     387           2 :     addComponent("Sxz");
     388           2 :     componentIsNotPeriodic("Sxz");
     389           2 :     addComponent("Syz");
     390           2 :     componentIsNotPeriodic("Syz");
     391             :   }
     392             : 
     393          29 :   requestAtoms(atoms, false);
     394          29 :   if(getDoScore()) {
     395          12 :     setParameters(coupl);
     396          12 :     Initialise(coupl.size());
     397             :   }
     398          29 :   setDerivatives();
     399          29 :   checkRead();
     400          29 : }
     401             : 
     402           5 : void RDC::do_svd() {
     403             : #ifdef __PLUMED_HAS_GSL
     404           5 :   gsl_vector_unique_ptr rdc_vec(gsl_vector_alloc(coupl.size())),
     405           5 :                         S(gsl_vector_alloc(5)),
     406           5 :                         Stmp(gsl_vector_alloc(5)),
     407           5 :                         work(gsl_vector_alloc(5)),
     408           5 :                         bc(gsl_vector_alloc(coupl.size()));
     409             : 
     410           5 :   gsl_matrix_unique_ptr coef_mat(gsl_matrix_alloc(coupl.size(),5)),
     411           5 :                         A(gsl_matrix_alloc(coupl.size(),5)),
     412           5 :                         V(gsl_matrix_alloc(5,5));
     413             : 
     414           5 :   gsl_matrix_set_zero(coef_mat.get());
     415           5 :   gsl_vector_set_zero(bc.get());
     416             : 
     417             :   unsigned index=0;
     418           5 :   std::vector<double> dmax(coupl.size());
     419          30 :   for(unsigned r=0; r<getNumberOfAtoms(); r+=2) {
     420          25 :     Vector  distance;
     421          25 :     if(pbc) {
     422          25 :       distance = pbcDistance(getPosition(r),getPosition(r+1));
     423             :     } else {
     424           0 :       distance = delta(getPosition(r),getPosition(r+1));
     425             :     }
     426          25 :     double d    = distance.modulo();
     427          25 :     double d2   = d*d;
     428          25 :     double d3   = d2*d;
     429          25 :     double id3  = 1./d3;
     430          25 :     double max  = -Const*mu_s*scale;
     431          25 :     dmax[index] = id3*max;
     432          25 :     double mu_x = distance[0]/d;
     433          25 :     double mu_y = distance[1]/d;
     434          25 :     double mu_z = distance[2]/d;
     435          25 :     gsl_vector_set(rdc_vec.get(),index,coupl[index]/dmax[index]);
     436          25 :     gsl_matrix_set(coef_mat.get(),index,0,gsl_matrix_get(coef_mat.get(),index,0)+(mu_x*mu_x-mu_z*mu_z));
     437          25 :     gsl_matrix_set(coef_mat.get(),index,1,gsl_matrix_get(coef_mat.get(),index,1)+(mu_y*mu_y-mu_z*mu_z));
     438          25 :     gsl_matrix_set(coef_mat.get(),index,2,gsl_matrix_get(coef_mat.get(),index,2)+(2.0*mu_x*mu_y));
     439          25 :     gsl_matrix_set(coef_mat.get(),index,3,gsl_matrix_get(coef_mat.get(),index,3)+(2.0*mu_x*mu_z));
     440          25 :     gsl_matrix_set(coef_mat.get(),index,4,gsl_matrix_get(coef_mat.get(),index,4)+(2.0*mu_y*mu_z));
     441          25 :     index++;
     442             :   }
     443           5 :   gsl_matrix_memcpy(A.get(),coef_mat.get());
     444           5 :   gsl_linalg_SV_decomp(A.get(), V.get(), Stmp.get(), work.get());
     445           5 :   gsl_linalg_SV_solve(A.get(), V.get(), Stmp.get(), rdc_vec.get(), S.get());
     446             :   /* tensor */
     447             :   Value* tensor;
     448           5 :   tensor=getPntrToComponent("Sxx");
     449           5 :   double Sxx = gsl_vector_get(S.get(),0);
     450             :   tensor->set(Sxx);
     451           5 :   tensor=getPntrToComponent("Syy");
     452           5 :   double Syy = gsl_vector_get(S.get(),1);
     453             :   tensor->set(Syy);
     454           5 :   tensor=getPntrToComponent("Szz");
     455           5 :   double Szz = -Sxx-Syy;
     456             :   tensor->set(Szz);
     457           5 :   tensor=getPntrToComponent("Sxy");
     458           5 :   double Sxy = gsl_vector_get(S.get(),2);
     459             :   tensor->set(Sxy);
     460           5 :   tensor=getPntrToComponent("Sxz");
     461           5 :   double Sxz = gsl_vector_get(S.get(),3);
     462             :   tensor->set(Sxz);
     463           5 :   tensor=getPntrToComponent("Syz");
     464           5 :   double Syz = gsl_vector_get(S.get(),4);
     465             :   tensor->set(Syz);
     466             : 
     467           5 :   gsl_blas_dgemv(CblasNoTrans, 1.0, coef_mat.get(), S.get(), 0., bc.get());
     468          30 :   for(index=0; index<coupl.size(); index++) {
     469          25 :     double rdc = gsl_vector_get(bc.get(),index)*dmax[index];
     470          25 :     Value* val=getPntrToComponent(index);
     471             :     val->set(rdc);
     472             :   }
     473             : #endif
     474           5 : }
     475             : 
     476        2172 : void RDC::calculate() {
     477        2172 :   if(svd) {
     478           5 :     do_svd();
     479           5 :     return;
     480             :   }
     481             : 
     482        2167 :   const double max  = -Const*scale*mu_s;
     483             :   const unsigned N=getNumberOfAtoms();
     484        2167 :   std::vector<Vector> dRDC(N/2, Vector{0.,0.,0.});
     485             : 
     486             :   /* RDC Calculations and forces */
     487        2167 :   #pragma omp parallel num_threads(OpenMP::getNumThreads())
     488             :   {
     489             :     #pragma omp for
     490             :     for(unsigned r=0; r<N; r+=2) {
     491             :       const unsigned index=r/2;
     492             :       Vector  distance;
     493             :       if(pbc) {
     494             :         distance   = pbcDistance(getPosition(r),getPosition(r+1));
     495             :       } else {
     496             :         distance   = delta(getPosition(r),getPosition(r+1));
     497             :       }
     498             :       const double d2    = distance.modulo2();
     499             :       const double ind   = 1./std::sqrt(d2);
     500             :       const double ind2  = 1./d2;
     501             :       const double ind3  = ind2*ind;
     502             :       const double x2    = distance[0]*distance[0]*ind2;
     503             :       const double y2    = distance[1]*distance[1]*ind2;
     504             :       const double z2    = distance[2]*distance[2]*ind2;
     505             :       const double dmax  = ind3*max;
     506             :       const double ddmax = dmax*ind2;
     507             : 
     508             :       const double rdc   = 0.5*dmax*(3.*z2-1.);
     509             :       const double prod_xy = (x2+y2-4.*z2);
     510             :       const double prod_z =  (3.*x2 + 3.*y2 - 2.*z2);
     511             : 
     512             :       dRDC[index] = -1.5*ddmax*distance;
     513             :       dRDC[index][0] *= prod_xy;
     514             :       dRDC[index][1] *= prod_xy;
     515             :       dRDC[index][2] *= prod_z;
     516             : 
     517             :       std::string num;
     518             :       Tools::convert(index,num);
     519             :       Value* val=getPntrToComponent("rdc-"+num);
     520             :       val->set(rdc);
     521             :       if(!getDoScore()) {
     522             :         setBoxDerivatives(val, Tensor(distance,dRDC[index]));
     523             :         setAtomsDerivatives(val, r,  dRDC[index]);
     524             :         setAtomsDerivatives(val, r+1, -dRDC[index]);
     525             :       } else {
     526             :         setCalcData(index, rdc);
     527             :       }
     528             :     }
     529             :   }
     530             : 
     531        2167 :   if(getDoScore()) {
     532             :     /* Metainference */
     533        1824 :     Tensor dervir;
     534        1824 :     double score = getScore();
     535        1824 :     setScore(score);
     536             : 
     537             :     /* calculate final derivatives */
     538        1824 :     Value* val=getPntrToComponent("score");
     539        9120 :     for(unsigned r=0; r<N; r+=2) {
     540        7296 :       const unsigned index=r/2;
     541        7296 :       Vector  distance;
     542        7296 :       if(pbc) {
     543        7296 :         distance   = pbcDistance(getPosition(r),getPosition(r+1));
     544             :       } else {
     545           0 :         distance   = delta(getPosition(r),getPosition(r+1));
     546             :       }
     547        7296 :       const Vector der = dRDC[index]*getMetaDer(index);
     548        7296 :       dervir += Tensor(distance, der);
     549        7296 :       setAtomsDerivatives(val, r,  der);
     550        7296 :       setAtomsDerivatives(val, r+1, -der);
     551             :     }
     552        1824 :     setBoxDerivatives(val, dervir);
     553             :   }
     554             : }
     555             : 
     556         327 : void RDC::update() {
     557             :   // write status file
     558         327 :   if(getWstride()>0&& (getStep()%getWstride()==0 || getCPT()) ) {
     559          29 :     writeStatus();
     560             :   }
     561         327 : }
     562             : 
     563             : }
     564             : }

Generated by: LCOV version 1.16