LCOV - code coverage report
Current view: top level - colvar - EEFSolv.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 230 244 94.3 %
Date: 2024-10-18 14:00:25 Functions: 7 8 87.5 %

          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             : 
      23             : /* This class was originally written by Thomas Loehr */
      24             : 
      25             : #include "Colvar.h"
      26             : #include "core/ActionRegister.h"
      27             : #include "core/ActionSet.h"
      28             : #include "core/PlumedMain.h"
      29             : #include "core/GenericMolInfo.h"
      30             : #include "tools/Communicator.h"
      31             : #include "tools/OpenMP.h"
      32             : #include <initializer_list>
      33             : 
      34             : #define INV_PI_SQRT_PI 0.179587122
      35             : #define KCAL_TO_KJ 4.184
      36             : #define ANG_TO_NM 0.1
      37             : #define ANG3_TO_NM3 0.001
      38             : 
      39             : namespace PLMD {
      40             : namespace colvar {
      41             : 
      42             : //+PLUMEDOC COLVAR EEFSOLV
      43             : /*
      44             : Calculates EEF1 solvation free energy for a group of atoms.
      45             : 
      46             : EEF1 is a solvent-accessible surface area based model, where the free energy of solvation is computed using a pairwise interaction term for non-hydrogen atoms:
      47             : \f[
      48             :     \Delta G^\mathrm{solv}_i = \Delta G^\mathrm{ref}_i - \sum_{j \neq i} f_i(r_{ij}) V_j
      49             : \f]
      50             : where \f$\Delta G^\mathrm{solv}_i\f$ is the free energy of solvation, \f$\Delta G^\mathrm{ref}_i\f$ is the reference solvation free energy, \f$V_j\f$ is the volume of atom \f$j\f$ and
      51             : \f[
      52             :     f_i(r) 4\pi r^2 = \frac{2}{\sqrt{\pi}} \frac{\Delta G^\mathrm{free}_i}{\lambda_i} \exp\left\{ - \frac{(r-R_i)^2}{\lambda^2_i}\right\}
      53             : \f]
      54             : where \f$\Delta G^\mathrm{free}_i\f$ is the solvation free energy of the isolated group, \f$\lambda_i\f$ is the correlation length equal to the width of the first solvation shell and \f$R_i\f$ is the van der Waals radius of atom \f$i\f$.
      55             : 
      56             : The output from this collective variable, the free energy of solvation, can be used with the \ref BIASVALUE keyword to provide implicit solvation to a system. All parameters are designed to be used with a modified CHARMM36 force field. It takes only non-hydrogen atoms as input, these can be conveniently specified using the \ref GROUP action with the NDX_GROUP parameter. To speed up the calculation, EEFSOLV internally uses a neighbor list with a cutoff dependent on the type of atom (maximum of 1.95 nm). This cutoff can be extended further by using the NL_BUFFER keyword.
      57             : 
      58             : \par Examples
      59             : 
      60             : \plumedfile
      61             : #SETTINGS MOLFILE=regtest/basic/rt77/peptide.pdb
      62             : MOLINFO MOLTYPE=protein STRUCTURE=peptide.pdb
      63             : WHOLEMOLECULES ENTITY0=1-111
      64             : 
      65             : # This allows us to select only non-hydrogen atoms
      66             : #SETTINGS AUXFILE=regtest/basic/rt77/index.ndx
      67             : protein-h: GROUP NDX_FILE=index.ndx NDX_GROUP=Protein-H
      68             : 
      69             : # We extend the cutoff by 0.1 nm and update the neighbor list every 40 steps
      70             : solv: EEFSOLV ATOMS=protein-h
      71             : 
      72             : # Here we actually add our calculated energy back to the potential
      73             : bias: BIASVALUE ARG=solv
      74             : 
      75             : PRINT ARG=solv FILE=SOLV
      76             : \endplumedfile
      77             : 
      78             : */
      79             : //+ENDPLUMEDOC
      80             : 
      81             : class EEFSolv : public Colvar {
      82             : private:
      83             :   bool pbc;
      84             :   bool serial;
      85             :   double delta_g_ref;
      86             :   double nl_buffer;
      87             :   unsigned nl_stride;
      88             :   unsigned nl_update;
      89             :   std::vector<std::vector<unsigned> > nl;
      90             :   std::vector<std::vector<bool> > nlexpo;
      91             :   std::vector<std::vector<double> > parameter;
      92             :   void setupConstants(const std::vector<AtomNumber> &atoms, std::vector<std::vector<double> > &parameter, bool tcorr);
      93             :   std::map<std::string, std::map<std::string, std::string> > setupTypeMap();
      94             :   std::map<std::string, std::vector<double> > setupValueMap();
      95             :   void update_neighb();
      96             : 
      97             : public:
      98             :   static void registerKeywords(Keywords& keys);
      99             :   explicit EEFSolv(const ActionOptions&);
     100             :   void calculate() override;
     101             : };
     102             : 
     103             : PLUMED_REGISTER_ACTION(EEFSolv,"EEFSOLV")
     104             : 
     105           7 : void EEFSolv::registerKeywords(Keywords& keys) {
     106           7 :   Colvar::registerKeywords(keys);
     107          14 :   keys.add("atoms", "ATOMS", "The atoms to be included in the calculation, e.g. the whole protein.");
     108          14 :   keys.add("compulsory", "NL_BUFFER", "0.1", "The buffer to the intrinsic cutoff used when calculating pairwise interactions.");
     109          14 :   keys.add("compulsory", "NL_STRIDE", "40", "The frequency with which the neighbor list is updated.");
     110          14 :   keys.addFlag("SERIAL",false,"Perform the calculation in serial - for debug purpose");
     111          14 :   keys.addFlag("TEMP_CORRECTION", false, "Correct free energy of solvation constants for temperatures different from 298.15 K");
     112           7 :   keys.setValueDescription("the EEF1 solvation free energy for the input atoms");
     113           7 : }
     114             : 
     115           5 : EEFSolv::EEFSolv(const ActionOptions&ao):
     116             :   PLUMED_COLVAR_INIT(ao),
     117           5 :   pbc(true),
     118           5 :   serial(false),
     119           5 :   delta_g_ref(0.),
     120           5 :   nl_buffer(0.1),
     121           5 :   nl_stride(40),
     122           5 :   nl_update(0)
     123             : {
     124             :   std::vector<AtomNumber> atoms;
     125          10 :   parseAtomList("ATOMS", atoms);
     126             :   const unsigned size = atoms.size();
     127           5 :   bool tcorr = false;
     128           5 :   parseFlag("TEMP_CORRECTION", tcorr);
     129           5 :   parse("NL_BUFFER", nl_buffer);
     130           5 :   parse("NL_STRIDE", nl_stride);
     131             : 
     132           5 :   bool nopbc = !pbc;
     133           5 :   parseFlag("NOPBC", nopbc);
     134           5 :   pbc = !nopbc;
     135             : 
     136           5 :   parseFlag("SERIAL",serial);
     137             : 
     138           5 :   checkRead();
     139             : 
     140          10 :   log << "  Bibliography " << plumed.cite("Lazaridis T, Karplus M, Proteins Struct. Funct. Genet. 35, 133 (1999)"); log << "\n";
     141             : 
     142           5 :   nl.resize(size);
     143           5 :   nlexpo.resize(size);
     144           5 :   parameter.resize(size, std::vector<double>(4, 0));
     145           5 :   setupConstants(atoms, parameter, tcorr);
     146             : 
     147           5 :   addValueWithDerivatives();
     148           5 :   setNotPeriodic();
     149           5 :   requestAtoms(atoms);
     150           5 : }
     151             : 
     152          30 : void EEFSolv::update_neighb() {
     153             :   const double lower_c2 = 0.24 * 0.24; // this is the cut-off for bonded atoms
     154             :   const unsigned size = getNumberOfAtoms();
     155             : 
     156        1830 :   for (unsigned i=0; i<size; i++) {
     157        1800 :     nl[i].clear();
     158             :     nlexpo[i].clear();
     159        1800 :     const Vector posi = getPosition(i);
     160             :     // Loop through neighboring atoms, add the ones below cutoff
     161       54900 :     for (unsigned j=i+1; j<size; j++) {
     162       53100 :       if(parameter[i][1]==0&&parameter[j][1]==0) continue;
     163       51750 :       const double d2 = delta(posi, getPosition(j)).modulo2();
     164       51750 :       if (d2 < lower_c2 && j < i+14) {
     165             :         // crude approximation for i-i+1/2 interactions,
     166             :         // we want to exclude atoms separated by less than three bonds
     167        2695 :         continue;
     168             :       }
     169             :       // We choose the maximum lambda value and use a more conservative cutoff
     170       49055 :       double mlambda = 1./parameter[i][2];
     171       49055 :       if (1./parameter[j][2] > mlambda) mlambda = 1./parameter[j][2];
     172       49055 :       const double c2 = (2. * mlambda + nl_buffer) * (2. * mlambda + nl_buffer);
     173       49055 :       if (d2 < c2 ) {
     174       26069 :         nl[i].push_back(j);
     175       26069 :         if(parameter[i][2] == parameter[j][2] && parameter[i][3] == parameter[j][3]) {
     176        5175 :           nlexpo[i].push_back(true);
     177       20894 :         } else nlexpo[i].push_back(false);
     178             :       }
     179             :     }
     180             :   }
     181          30 : }
     182             : 
     183          30 : void EEFSolv::calculate() {
     184          30 :   if(pbc) makeWhole();
     185          30 :   if(getExchangeStep()) nl_update = 0;
     186          30 :   if(nl_update==0) update_neighb();
     187             : 
     188             :   const unsigned size=getNumberOfAtoms();
     189          30 :   double bias = 0.0;
     190          30 :   std::vector<Vector> deriv(size, Vector(0,0,0));
     191             : 
     192             :   unsigned stride;
     193             :   unsigned rank;
     194          30 :   if(serial) {
     195             :     stride=1;
     196             :     rank=0;
     197             :   } else {
     198          30 :     stride=comm.Get_size();
     199          30 :     rank=comm.Get_rank();
     200             :   }
     201             : 
     202          30 :   unsigned nt=OpenMP::getNumThreads();
     203          30 :   if(nt*stride*10>size) nt=1;
     204             : 
     205          30 :   #pragma omp parallel num_threads(nt)
     206             :   {
     207             :     std::vector<Vector> deriv_omp(size, Vector(0,0,0));
     208             :     #pragma omp for reduction(+:bias) nowait
     209             :     for (unsigned i=rank; i<size; i+=stride) {
     210             :       const Vector posi = getPosition(i);
     211             :       double fedensity = 0.0;
     212             :       Vector deriv_i;
     213             :       const double vdw_volume_i   = parameter[i][0];
     214             :       const double delta_g_free_i = parameter[i][1];
     215             :       const double inv_lambda_i   = parameter[i][2];
     216             :       const double vdw_radius_i   = parameter[i][3];
     217             : 
     218             :       // The pairwise interactions are unsymmetric, but we can get away with calculating the distance only once
     219             :       for (unsigned i_nl=0; i_nl<nl[i].size(); i_nl++) {
     220             :         const unsigned j = nl[i][i_nl];
     221             :         const double vdw_volume_j   = parameter[j][0];
     222             :         const double delta_g_free_j = parameter[j][1];
     223             :         const double inv_lambda_j   = parameter[j][2];
     224             :         const double vdw_radius_j   = parameter[j][3];
     225             : 
     226             :         const Vector dist     = delta(posi, getPosition(j));
     227             :         const double rij      = dist.modulo();
     228             :         const double inv_rij  = 1.0 / rij;
     229             :         const double inv_rij2 = inv_rij * inv_rij;
     230             :         const double fact_ij  = inv_rij2 * delta_g_free_i * vdw_volume_j * INV_PI_SQRT_PI * inv_lambda_i;
     231             :         const double fact_ji  = inv_rij2 * delta_g_free_j * vdw_volume_i * INV_PI_SQRT_PI * inv_lambda_j;
     232             : 
     233             :         // in this case we can calculate a single exponential
     234             :         if(!nlexpo[i][i_nl]) {
     235             :           // i-j interaction
     236             :           if(inv_rij > 0.5*inv_lambda_i && delta_g_free_i!=0.)
     237             :           {
     238             :             const double e_arg = (rij - vdw_radius_i)*inv_lambda_i;
     239             :             const double expo  = std::exp(-e_arg*e_arg);
     240             :             const double fact  = expo*fact_ij;
     241             :             const double e_deriv = inv_rij*fact*(inv_rij + e_arg*inv_lambda_i);
     242             :             const Vector dd    = e_deriv*dist;
     243             :             fedensity    += fact;
     244             :             deriv_i      += dd;
     245             :             if(nt>1) deriv_omp[j] -= dd;
     246             :             else deriv[j] -= dd;
     247             :           }
     248             : 
     249             :           // j-i interaction
     250             :           if(inv_rij > 0.5*inv_lambda_j && delta_g_free_j!=0.)
     251             :           {
     252             :             const double e_arg = (rij - vdw_radius_j)*inv_lambda_j;
     253             :             const double expo  = std::exp(-e_arg*e_arg);
     254             :             const double fact  = expo*fact_ji;
     255             :             const double e_deriv = inv_rij*fact*(inv_rij + e_arg*inv_lambda_j);
     256             :             const Vector dd    = e_deriv*dist;
     257             :             fedensity    += fact;
     258             :             deriv_i      += dd;
     259             :             if(nt>1) deriv_omp[j] -= dd;
     260             :             else deriv[j] -= dd;
     261             :           }
     262             :         } else {
     263             :           // i-j interaction
     264             :           if(inv_rij > 0.5*inv_lambda_i)
     265             :           {
     266             :             const double e_arg = (rij - vdw_radius_i)*inv_lambda_i;
     267             :             const double expo  = std::exp(-e_arg*e_arg);
     268             :             const double fact  = expo*(fact_ij + fact_ji);
     269             :             const double e_deriv = inv_rij*fact*(inv_rij + e_arg*inv_lambda_i);
     270             :             const Vector dd    = e_deriv*dist;
     271             :             fedensity    += fact;
     272             :             deriv_i      += dd;
     273             :             if(nt>1) deriv_omp[j] -= dd;
     274             :             else deriv[j] -= dd;
     275             :           }
     276             :         }
     277             : 
     278             :       }
     279             :       if(nt>1) deriv_omp[i] += deriv_i;
     280             :       else deriv[i] += deriv_i;
     281             :       bias += 0.5*fedensity;
     282             :     }
     283             :     #pragma omp critical
     284             :     if(nt>1) for(unsigned i=0; i<size; i++) deriv[i]+=deriv_omp[i];
     285             :   }
     286             : 
     287          30 :   if(!serial) {
     288          30 :     comm.Sum(bias);
     289          30 :     if(!deriv.empty()) comm.Sum(&deriv[0][0],3*deriv.size());
     290             :   }
     291             : 
     292          30 :   Tensor virial;
     293        1830 :   for(unsigned i=0; i<size; i++) {
     294        1800 :     setAtomsDerivatives(i, -deriv[i]);
     295        1800 :     virial += Tensor(getPosition(i), -deriv[i]);
     296             :   }
     297          30 :   setBoxDerivatives(-virial);
     298          30 :   setValue(delta_g_ref - bias);
     299             : 
     300             :   // Keep track of the neighbourlist updates
     301          30 :   nl_update++;
     302          30 :   if (nl_update == nl_stride) {
     303          30 :     nl_update = 0;
     304             :   }
     305          30 : }
     306             : 
     307           5 : void EEFSolv::setupConstants(const std::vector<AtomNumber> &atoms, std::vector<std::vector<double> > &parameter, bool tcorr) {
     308             :   std::vector<std::vector<double> > parameter_temp;
     309          10 :   parameter_temp.resize(atoms.size(), std::vector<double>(7,0));
     310             :   std::map<std::string, std::vector<double> > valuemap;
     311             :   std::map<std::string, std::map<std::string, std::string> > typemap;
     312           5 :   valuemap = setupValueMap();
     313           5 :   typemap  = setupTypeMap();
     314           5 :   auto * moldat = plumed.getActionSet().selectLatest<GenericMolInfo*>(this);
     315             :   bool cter=false;
     316           5 :   if (moldat) {
     317           5 :     log<<"  MOLINFO DATA found with label " <<moldat->getLabel()<<", using proper atom names\n";
     318         305 :     for(unsigned i=0; i<atoms.size(); ++i) {
     319             : 
     320             :       // Get atom and residue names
     321         300 :       std::string Aname = moldat->getAtomName(atoms[i]);
     322         300 :       std::string Rname = moldat->getResidueName(atoms[i]);
     323         300 :       std::string Atype = typemap[Rname][Aname];
     324             : 
     325             :       // Check for terminal COOH or COO- (different atomtypes & parameters!)
     326         598 :       if (Aname == "OT1" || Aname == "OXT") {
     327             :         // We create a temporary AtomNumber object to access future atoms
     328             :         unsigned ai = atoms[i].index();
     329             :         AtomNumber tmp_an;
     330           2 :         tmp_an.setIndex(ai + 2);
     331           2 :         if (moldat->checkForAtom(tmp_an) && moldat->getAtomName(tmp_an) == "HT2") {
     332             :           // COOH
     333             :           Atype = "OB";
     334             :         } else {
     335             :           // COO-
     336             :           Atype = "OC";
     337             :         }
     338             :         cter = true;
     339             :       }
     340         302 :       if (Aname == "OT2" || (cter == true && Aname == "O")) {
     341             :         unsigned ai = atoms[i].index();
     342             :         AtomNumber tmp_an;
     343           2 :         tmp_an.setIndex(ai + 1);
     344           2 :         if (moldat->checkForAtom(tmp_an) && moldat->getAtomName(tmp_an) == "HT2") {
     345             :           // COOH
     346             :           Atype = "OH1";
     347             :         } else {
     348             :           // COO-
     349             :           Atype = "OC";
     350             :         }
     351             :       }
     352             : 
     353             :       // Check for H-atoms
     354             :       char type;
     355         300 :       char first = Aname.at(0);
     356             : 
     357             :       // GOLDEN RULE: type is first letter, if not a number
     358         300 :       if (!isdigit(first)) {
     359             :         type = first;
     360             :         // otherwise is the second
     361             :       } else {
     362           0 :         type = Aname.at(1);
     363             :       }
     364             : 
     365         300 :       if (type == 'H') {
     366           0 :         error("EEF1-SB does not allow the use of hydrogen atoms!\n");
     367             :       }
     368             : 
     369             :       // Lookup atomtype in table or throw exception if its not there
     370             :       try {
     371         300 :         parameter_temp[i] = valuemap.at(Atype);
     372           0 :       } catch (const std::exception &e) {
     373           0 :         log << "Type: " << Atype << "  Name: " << Aname << "  Residue: " << Rname << "\n";
     374           0 :         error("Invalid atom type!\n");
     375           0 :       }
     376             : 
     377             :       // Temperature correction
     378         300 :       if (tcorr && parameter[i][1] > 0.0) {
     379             :         const double t0 = 298.15;
     380           0 :         const double delta_g_ref_t0 = parameter_temp[i][1];
     381           0 :         const double delta_h_ref_t0 = parameter_temp[i][3];
     382           0 :         const double delta_cp = parameter_temp[i][4];
     383           0 :         const double delta_s_ref_t0 = (delta_h_ref_t0 - delta_g_ref_t0) / t0;
     384           0 :         const double t = getkBT() / getKBoltzmann();
     385           0 :         parameter_temp[i][1] -= delta_s_ref_t0 * (t - t0) - delta_cp * t * std::log(t / t0) + delta_cp * (t - t0);
     386           0 :         parameter_temp[i][2] *= parameter_temp[i][1] / delta_g_ref_t0;
     387             :       }
     388         300 :       parameter[i][0] = parameter_temp[i][0];
     389         300 :       parameter[i][1] = parameter_temp[i][2];
     390         300 :       parameter[i][2] = parameter_temp[i][5];
     391         300 :       parameter[i][3] = parameter_temp[i][6];
     392             :     }
     393             :   } else {
     394           0 :     error("MOLINFO DATA not found\n");
     395             :   }
     396         305 :   for(unsigned i=0; i<atoms.size(); ++i) delta_g_ref += parameter_temp[i][1];
     397           5 : }
     398             : 
     399           5 : std::map<std::string, std::map<std::string, std::string> > EEFSolv::setupTypeMap()  {
     400             :   std::map<std::string, std::map<std::string, std::string> > typemap;
     401          10 :   typemap["ACE"] = {
     402             :     {"CH3", "CT3"},
     403             :     {"HH31","HA3"},
     404             :     {"HH32","HA3"},
     405             :     {"HH33","HA3"},
     406             :     {"C",   "C"  },
     407             :     {"O",   "O"  }
     408          40 :   };
     409          10 :   typemap["ALA"] = {
     410             :     {"N",   "NH1"},
     411             :     {"HN",  "H"  },
     412             :     {"CA",  "CT1"},
     413             :     {"HA",  "HB1"},
     414             :     {"CB",  "CT3"},
     415             :     {"HB1", "HA3"},
     416             :     {"HB2", "HA3"},
     417             :     {"HB3", "HA3"},
     418             :     {"C",   "C"  },
     419             :     {"O",   "O"  }
     420          60 :   };
     421          10 :   typemap["ARG"] = {
     422             :     {"N",    "NH1"},
     423             :     {"HN",   "H"  },
     424             :     {"CA",   "CT1"},
     425             :     {"HA",   "HB1"},
     426             :     {"CB",   "CT2"},
     427             :     {"HB1",  "HA2"},
     428             :     {"HB2",  "HA2"},
     429             :     {"CG",   "CT2"},
     430             :     {"HG1",  "HA2"},
     431             :     {"HG2",  "HA2"},
     432             :     {"CD",   "CT2"},
     433             :     {"HD1",  "HA2"},
     434             :     {"HD2",  "HA2"},
     435             :     {"NE",   "NC2"},
     436             :     {"HE",   "HC" },
     437             :     {"CZ",   "C"  },
     438             :     {"NH1",  "NC2"},
     439             :     {"HH11", "HC" },
     440             :     {"HH12", "HC" },
     441             :     {"NH2",  "NC2"},
     442             :     {"HH21", "HC" },
     443             :     {"HH22", "HC" },
     444             :     {"C",    "C"  },
     445             :     {"O",    "O"  }
     446         130 :   };
     447          10 :   typemap["ASN"] = {
     448             :     {"N",    "NH1"},
     449             :     {"HN",   "H"  },
     450             :     {"CA",   "CT1"},
     451             :     {"HA",   "HB1"},
     452             :     {"CB",   "CT2"},
     453             :     {"HB1",  "HA2"},
     454             :     {"HB2",  "HA2"},
     455             :     {"CG",   "CC" },
     456             :     {"OD1",  "O"  },
     457             :     {"ND2",  "NH2"},
     458             :     {"HD21", "H"  },
     459             :     {"HD22", "H"  },
     460             :     {"C",    "C"  },
     461             :     {"O",    "O"  }
     462          80 :   };
     463          10 :   typemap["ASPP"] = {
     464             :     {"N",   "NH1"},
     465             :     {"HN",  "H"  },
     466             :     {"CA",  "CT1"},
     467             :     {"HA",  "HB1"},
     468             :     {"CB",  "CT2"},
     469             :     {"HB1", "HA2"},
     470             :     {"HB2", "HA2"},
     471             :     {"CG",  "CD" },
     472             :     {"OD1", "OB" },
     473             :     {"OD2", "OH1"},
     474             :     {"HD2", "H"  },
     475             :     {"C",   "C"  },
     476             :     {"O",   "O"  }
     477          75 :   };
     478          10 :   typemap["ASP"] = {
     479             :     {"N",   "NH1"},
     480             :     {"HN",  "H"  },
     481             :     {"CA",  "CT1"},
     482             :     {"HA",  "HB1"},
     483             :     {"CB",  "CT2"},
     484             :     {"HB1", "HA2"},
     485             :     {"HB2", "HA2"},
     486             :     {"CG",  "CC" },
     487             :     {"OD1", "OC" },
     488             :     {"OD2", "OC" },
     489             :     {"C",   "C"  },
     490             :     {"O",   "O"  }
     491          70 :   };
     492          10 :   typemap["CYS"] = {
     493             :     {"N",   "NH1"},
     494             :     {"HN",  "H"  },
     495             :     {"CA",  "CT1"},
     496             :     {"HA",  "HB1"},
     497             :     {"CB",  "CT2"},
     498             :     {"HB1", "HA2"},
     499             :     {"HB2", "HA2"},
     500             :     {"SG",  "S"  },
     501             :     {"HG1", "HS" },
     502             :     {"C",   "C"  },
     503             :     {"O",   "O"  }
     504          65 :   };
     505          10 :   typemap["GLN"] = {
     506             :     {"N",    "NH1" },
     507             :     {"HN",   "H"   },
     508             :     {"CA",   "CT1" },
     509             :     {"HA",   "HB1" },
     510             :     {"CB",   "CT2" },
     511             :     {"HB1",  "HA2" },
     512             :     {"HB2",  "HA2" },
     513             :     {"CG",   "CT2" },
     514             :     {"HG1",  "HA2" },
     515             :     {"HG2",  "HA2" },
     516             :     {"CD",   "CC"  },
     517             :     {"OE1",  "O"   },
     518             :     {"NE2",  "NH2" },
     519             :     {"HE21", "H"   },
     520             :     {"HE22", "H"   },
     521             :     {"C",    "C"   },
     522             :     {"O",    "O"   }
     523          95 :   };
     524          10 :   typemap["GLUP"] = {
     525             :     {"N",   "NH1"},
     526             :     {"HN",  "H"  },
     527             :     {"CA",  "CT1"},
     528             :     {"HA",  "HB1"},
     529             :     {"CB",  "CT2"},
     530             :     {"HB1", "HA2"},
     531             :     {"HB2", "HA2"},
     532             :     {"CG",  "CT2"},
     533             :     {"HG1", "HA2"},
     534             :     {"HG2", "HA2"},
     535             :     {"CD",  "CD" },
     536             :     {"OE1", "OB" },
     537             :     {"OE2", "OH1"},
     538             :     {"HE2", "H"  },
     539             :     {"C",   "C"  },
     540             :     {"O",   "O"  }
     541          90 :   };
     542          10 :   typemap["GLU"] = {
     543             :     {"N",   "NH1"},
     544             :     {"HN",  "H"  },
     545             :     {"CA",  "CT1"},
     546             :     {"HA",  "HB1"},
     547             :     {"CB",  "CT2"},
     548             :     {"HB1", "HA2"},
     549             :     {"HB2", "HA2"},
     550             :     {"CG",  "CT2"},
     551             :     {"HG1", "HA2"},
     552             :     {"HG2", "HA2"},
     553             :     {"CD",  "CC" },
     554             :     {"OE1", "OC" },
     555             :     {"OE2", "OC" },
     556             :     {"C",   "C"  },
     557             :     {"O",   "O"  }
     558          85 :   };
     559          10 :   typemap["GLY"] = {
     560             :     {"N",   "NH1"},
     561             :     {"HN",  "H"  },
     562             :     {"CA",  "CT2"},
     563             :     {"HA1", "HB2"},
     564             :     {"HA2", "HB2"},
     565             :     {"C",   "C"  },
     566             :     {"O",   "O"  }
     567          45 :   };
     568          10 :   typemap["HSD"] = {
     569             :     {"N",   "NH1"},
     570             :     {"HN",  "H"  },
     571             :     {"CA",  "CT1"},
     572             :     {"HA",  "HB1"},
     573             :     {"CB",  "CT2"},
     574             :     {"HB1", "HA2"},
     575             :     {"HB2", "HA2"},
     576             :     {"ND1", "NR1"},
     577             :     {"HD1", "H"  },
     578             :     {"CG",  "CPH1"},
     579             :     {"CE1", "CPH2"},
     580             :     {"HE1", "HR1"},
     581             :     {"NE2", "NR2"},
     582             :     {"CD2", "CPH1"},
     583             :     {"HD2", "HR3"},
     584             :     {"C",   "C"  },
     585             :     {"O",   "O"  }
     586          95 :   };
     587          10 :   typemap["HIS"] = {
     588             :     {"N",   "NH1"},
     589             :     {"HN",  "H"  },
     590             :     {"CA",  "CT1"},
     591             :     {"HA",  "HB1"},
     592             :     {"CB",  "CT2"},
     593             :     {"HB1", "HA2"},
     594             :     {"HB2", "HA2"},
     595             :     {"ND1", "NR2"},
     596             :     {"CG",  "CPH1"},
     597             :     {"CE1", "CPH2"},
     598             :     {"HE1", "HR1"},
     599             :     {"NE2", "NR1"},
     600             :     {"HE2", "H"  },
     601             :     {"CD2", "CPH1"},
     602             :     {"HD2", "HR3"},
     603             :     {"C",   "C"  },
     604             :     {"O",   "O"  }
     605          95 :   };
     606          10 :   typemap["HSE"] = {
     607             :     {"N",   "NH1"},
     608             :     {"HN",  "H"  },
     609             :     {"CA",  "CT1"},
     610             :     {"HA",  "HB1"},
     611             :     {"CB",  "CT2"},
     612             :     {"HB1", "HA2"},
     613             :     {"HB2", "HA2"},
     614             :     {"ND1", "NR2"},
     615             :     {"CG",  "CPH1"},
     616             :     {"CE1", "CPH2"},
     617             :     {"HE1", "HR1"},
     618             :     {"NE2", "NR1"},
     619             :     {"HE2", "H"  },
     620             :     {"CD2", "CPH1"},
     621             :     {"HD2", "HR3"},
     622             :     {"C",   "C"  },
     623             :     {"O",   "O"  }
     624          95 :   };
     625          10 :   typemap["HSP"] = {
     626             :     {"N",   "NH1"},
     627             :     {"HN",  "H"  },
     628             :     {"CA",  "CT1"},
     629             :     {"HA",  "HB1"},
     630             :     {"CB",  "CT2"},
     631             :     {"HB1", "HA2"},
     632             :     {"HB2", "HA2"},
     633             :     {"CD2", "CPH1"},
     634             :     {"HD2", "HR1"},
     635             :     {"CG",  "CPH1"},
     636             :     {"NE2", "NR3"},
     637             :     {"HE2", "H"  },
     638             :     {"ND1", "NR3"},
     639             :     {"HD1", "H"  },
     640             :     {"CE1", "CPH2"},
     641             :     {"HE1", "HR2"},
     642             :     {"C",   "C"  },
     643             :     {"O",   "O"  }
     644         100 :   };
     645          10 :   typemap["ILE"] = {
     646             :     {"N",    "NH1"},
     647             :     {"HN",   "H"  },
     648             :     {"CA",   "CT1"},
     649             :     {"HA",   "HB1"},
     650             :     {"CB",   "CT1"},
     651             :     {"HB",   "HA1"},
     652             :     {"CG2",  "CT3"},
     653             :     {"HG21", "HA3"},
     654             :     {"HG22", "HA3"},
     655             :     {"HG23", "HA3"},
     656             :     {"CG1",  "CT2"},
     657             :     {"HG11", "HA2"},
     658             :     {"HG12", "HA2"},
     659             :     {"CD",   "CT3"},
     660             :     {"HD1",  "HA3"},
     661             :     {"HD2",  "HA3"},
     662             :     {"HD3",  "HA3"},
     663             :     {"C",    "C"  },
     664             :     {"O",    "O"  }
     665         105 :   };
     666          10 :   typemap["LEU"] = {
     667             :     {"N",    "NH1"},
     668             :     {"HN",   "H"  },
     669             :     {"CA",   "CT1"},
     670             :     {"HA",   "HB1"},
     671             :     {"CB",   "CT2"},
     672             :     {"HB1",  "HA2"},
     673             :     {"HB2",  "HA2"},
     674             :     {"CG",   "CT1"},
     675             :     {"HG",   "HA1"},
     676             :     {"CD1",  "CT3"},
     677             :     {"HD11", "HA3"},
     678             :     {"HD12", "HA3"},
     679             :     {"HD13", "HA3"},
     680             :     {"CD2",  "CT3"},
     681             :     {"HD21", "HA3"},
     682             :     {"HD22", "HA3"},
     683             :     {"HD23", "HA3"},
     684             :     {"C",    "C"  },
     685             :     {"O",    "O"  }
     686         105 :   };
     687          10 :   typemap["LYS"] = {
     688             :     {"N",   "NH1"},
     689             :     {"HN",  "H"  },
     690             :     {"CA",  "CT1"},
     691             :     {"HA",  "HB1"},
     692             :     {"CB",  "CT2"},
     693             :     {"HB1", "HA2"},
     694             :     {"HB2", "HA2"},
     695             :     {"CG",  "CT2"},
     696             :     {"HG1", "HA2"},
     697             :     {"HG2", "HA2"},
     698             :     {"CD",  "CT2"},
     699             :     {"HD1", "HA2"},
     700             :     {"HD2", "HA2"},
     701             :     {"CE",  "CT2"},
     702             :     {"HE1", "HA2"},
     703             :     {"HE2", "HA2"},
     704             :     {"NZ",  "NH3"},
     705             :     {"HZ1", "HC" },
     706             :     {"HZ2", "HC" },
     707             :     {"HZ3", "HC" },
     708             :     {"C",   "C"  },
     709             :     {"O",   "O"  }
     710         120 :   };
     711          10 :   typemap["MET"] = {
     712             :     {"N",   "NH1"},
     713             :     {"HN",  "H"  },
     714             :     {"CA",  "CT1"},
     715             :     {"HA",  "HB1"},
     716             :     {"CB",  "CT2"},
     717             :     {"HB1", "HA2"},
     718             :     {"HB2", "HA2"},
     719             :     {"CG",  "CT2"},
     720             :     {"HG1", "HA2"},
     721             :     {"HG2", "HA2"},
     722             :     {"SD",  "S"  },
     723             :     {"CE",  "CT3"},
     724             :     {"HE1", "HA3"},
     725             :     {"HE2", "HA3"},
     726             :     {"HE3", "HA3"},
     727             :     {"C",   "C"  },
     728             :     {"O",   "O"  }
     729          95 :   };
     730          10 :   typemap["NMA"] = {
     731             :     {"N",   "NH1"},
     732             :     {"HN",  "H"  },
     733             :     {"CH3", "CT3"},
     734             :     {"HH31","HA3"},
     735             :     {"HH32","HA3"},
     736             :     {"HH33","HA3"},
     737          40 :   };
     738          10 :   typemap["PHE"] = {
     739             :     {"N",   "NH1"},
     740             :     {"HN",  "H"  },
     741             :     {"CA",  "CT1"},
     742             :     {"HA",  "HB1"},
     743             :     {"CB",  "CT2"},
     744             :     {"HB1", "HA2"},
     745             :     {"HB2", "HA2"},
     746             :     {"CG",  "CA" },
     747             :     {"CD1", "CA" },
     748             :     {"HD1", "HP" },
     749             :     {"CE1", "CA" },
     750             :     {"HE1", "HP" },
     751             :     {"CZ",  "CA" },
     752             :     {"HZ",  "HP" },
     753             :     {"CD2", "CA" },
     754             :     {"HD2", "HP" },
     755             :     {"CE2", "CA" },
     756             :     {"HE2", "HP" },
     757             :     {"C",   "C"  },
     758             :     {"O",   "O"  }
     759         110 :   };
     760          10 :   typemap["PRO"] = {
     761             :     {"N",   "N"  },
     762             :     {"CD",  "CP3"},
     763             :     {"HD1", "HA2"},
     764             :     {"HD2", "HA2"},
     765             :     {"CA",  "CP1"},
     766             :     {"HA",  "HB1"},
     767             :     {"CB",  "CP2"},
     768             :     {"HB1", "HA2"},
     769             :     {"HB2", "HA2"},
     770             :     {"CG",  "CP2"},
     771             :     {"HG1", "HA2"},
     772             :     {"HG2", "HA2"},
     773             :     {"C",   "C"  },
     774             :     {"O",   "O"  }
     775          80 :   };
     776          10 :   typemap["SER"] = {
     777             :     {"N",   "NH1"},
     778             :     {"HN",  "H"  },
     779             :     {"CA",  "CT1"},
     780             :     {"HA",  "HB1"},
     781             :     {"CB",  "CT2"},
     782             :     {"HB1", "HA2"},
     783             :     {"HB2", "HA2"},
     784             :     {"OG",  "OH1"},
     785             :     {"HG1", "H"  },
     786             :     {"C",   "C"  },
     787             :     {"O",   "O"  }
     788          65 :   };
     789          10 :   typemap["THR"] = {
     790             :     {"N",    "NH1"},
     791             :     {"HN",   "H"  },
     792             :     {"CA",   "CT1"},
     793             :     {"HA",   "HB1"},
     794             :     {"CB",   "CT1"},
     795             :     {"HB",   "HA1"},
     796             :     {"OG1",  "OH1"},
     797             :     {"HG1",  "H"  },
     798             :     {"CG2",  "CT3"},
     799             :     {"HG21", "HA3"},
     800             :     {"HG22", "HA3"},
     801             :     {"HG23", "HA3"},
     802             :     {"C",    "C"  },
     803             :     {"O",    "O"  }
     804          80 :   };
     805          10 :   typemap["TRP"] = {
     806             :     {"N",   "NH1"},
     807             :     {"HN",  "H"  },
     808             :     {"CA",  "CT1"},
     809             :     {"HA",  "HB1"},
     810             :     {"CB",  "CT2"},
     811             :     {"HB1", "HA2"},
     812             :     {"HB2", "HA2"},
     813             :     {"CG",  "CY" },
     814             :     {"CD1", "CA" },
     815             :     {"HD1", "HP" },
     816             :     {"NE1", "NY" },
     817             :     {"HE1", "H"  },
     818             :     {"CE2", "CPT"},
     819             :     {"CD2", "CPT"},
     820             :     {"CE3", "CAI"},
     821             :     {"HE3", "HP" },
     822             :     {"CZ3", "CA" },
     823             :     {"HZ3", "HP" },
     824             :     {"CZ2", "CAI"},
     825             :     {"HZ2", "HP" },
     826             :     {"CH2", "CA" },
     827             :     {"HH2", "HP" },
     828             :     {"C",   "C"  },
     829             :     {"O",   "O"  }
     830         130 :   };
     831          10 :   typemap["TYR"] = {
     832             :     {"N",   "NH1"},
     833             :     {"HN",  "H"  },
     834             :     {"CA",  "CT1"},
     835             :     {"HA",  "HB1"},
     836             :     {"CB",  "CT2"},
     837             :     {"HB1", "HA2"},
     838             :     {"HB2", "HA2"},
     839             :     {"CG",  "CA" },
     840             :     {"CD1", "CA" },
     841             :     {"HD1", "HP" },
     842             :     {"CE1", "CA" },
     843             :     {"HE1", "HP" },
     844             :     {"CZ",  "CA" },
     845             :     {"OH",  "OH1"},
     846             :     {"HH",  "H"  },
     847             :     {"CD2", "CA" },
     848             :     {"HD2", "HP" },
     849             :     {"CE2", "CA" },
     850             :     {"HE2", "HP" },
     851             :     {"C",   "C"  },
     852             :     {"O",   "O"  }
     853         115 :   };
     854          10 :   typemap["VAL"] = {
     855             :     {"N",    "NH1"},
     856             :     {"HN",   "H"  },
     857             :     {"CA",   "CT1"},
     858             :     {"HA",   "HB1"},
     859             :     {"CB",   "CT1"},
     860             :     {"HB",   "HA1"},
     861             :     {"CG1",  "CT3"},
     862             :     {"HG11", "HA3"},
     863             :     {"HG12", "HA3"},
     864             :     {"HG13", "HA3"},
     865             :     {"CG2",  "CT3"},
     866             :     {"HG21", "HA3"},
     867             :     {"HG22", "HA3"},
     868             :     {"HG23", "HA3"},
     869             :     {"C",    "C"  },
     870             :     {"O",    "O"  }
     871          90 :   };
     872           5 :   return typemap;
     873             : }
     874             : 
     875           5 : std::map<std::string, std::vector<double> > EEFSolv::setupValueMap() {
     876             :   // Volume ∆Gref ∆Gfree ∆H ∆Cp λ vdw_radius
     877             :   std::map<std::string, std::vector<double> > valuemap;
     878           5 :   valuemap["C"] = {
     879             :     ANG3_TO_NM3 * 14.720,
     880             :     KCAL_TO_KJ * 0.000,
     881             :     KCAL_TO_KJ * 0.000,
     882             :     KCAL_TO_KJ * 0.000,
     883             :     KCAL_TO_KJ * 0.0,
     884             :     1. / (ANG_TO_NM * 3.5),
     885             :     0.20,
     886          10 :   };
     887           5 :   valuemap["CD"] = {
     888             :     ANG3_TO_NM3 * 14.720,
     889             :     KCAL_TO_KJ * 0.000,
     890             :     KCAL_TO_KJ * 0.000,
     891             :     KCAL_TO_KJ * 0.000,
     892             :     KCAL_TO_KJ * 0.0,
     893             :     1. / (ANG_TO_NM * 3.5),
     894             :     0.20,
     895          10 :   };
     896           5 :   valuemap["CT1"] = {
     897             :     ANG3_TO_NM3 * 11.507,
     898             :     KCAL_TO_KJ * -0.187,
     899             :     KCAL_TO_KJ * -0.187,
     900             :     KCAL_TO_KJ * 0.876,
     901             :     KCAL_TO_KJ * 0.0,
     902             :     1. / (ANG_TO_NM * 3.5),
     903             :     0.20,
     904          10 :   };
     905           5 :   valuemap["CT2"] = {
     906             :     ANG3_TO_NM3 * 18.850,
     907             :     KCAL_TO_KJ * 0.372,
     908             :     KCAL_TO_KJ * 0.372,
     909             :     KCAL_TO_KJ * -0.610,
     910             :     KCAL_TO_KJ * 18.6,
     911             :     1. / (ANG_TO_NM * 3.5),
     912             :     0.20,
     913          10 :   };
     914           5 :   valuemap["CT2A"] = {
     915             :     ANG3_TO_NM3 * 18.666,
     916             :     KCAL_TO_KJ * 0.372,
     917             :     KCAL_TO_KJ * 0.372,
     918             :     KCAL_TO_KJ * -0.610,
     919             :     KCAL_TO_KJ * 18.6,
     920             :     1. / (ANG_TO_NM * 3.5),
     921             :     0.20,
     922          10 :   };
     923           5 :   valuemap["CT3"] = {
     924             :     ANG3_TO_NM3 * 27.941,
     925             :     KCAL_TO_KJ * 1.089,
     926             :     KCAL_TO_KJ * 1.089,
     927             :     KCAL_TO_KJ * -1.779,
     928             :     KCAL_TO_KJ * 35.6,
     929             :     1. / (ANG_TO_NM * 3.5),
     930             :     0.204,
     931          10 :   };
     932           5 :   valuemap["CPH1"] = {
     933             :     ANG3_TO_NM3 * 5.275,
     934             :     KCAL_TO_KJ * 0.057,
     935             :     KCAL_TO_KJ * 0.080,
     936             :     KCAL_TO_KJ * -0.973,
     937             :     KCAL_TO_KJ * 6.9,
     938             :     1. / (ANG_TO_NM * 3.5),
     939             :     0.18,
     940          10 :   };
     941           5 :   valuemap["CPH2"] = {
     942             :     ANG3_TO_NM3 * 11.796,
     943             :     KCAL_TO_KJ * 0.057,
     944             :     KCAL_TO_KJ * 0.080,
     945             :     KCAL_TO_KJ * -0.973,
     946             :     KCAL_TO_KJ * 6.9,
     947             :     1. / (ANG_TO_NM * 3.5),
     948             :     0.18,
     949          10 :   };
     950           5 :   valuemap["CPT"] = {
     951             :     ANG3_TO_NM3 * 4.669,
     952             :     KCAL_TO_KJ * -0.890,
     953             :     KCAL_TO_KJ * -0.890,
     954             :     KCAL_TO_KJ * 2.220,
     955             :     KCAL_TO_KJ * 6.9,
     956             :     1. / (ANG_TO_NM * 3.5),
     957             :     0.186,
     958          10 :   };
     959           5 :   valuemap["CY"] = {
     960             :     ANG3_TO_NM3 * 10.507,
     961             :     KCAL_TO_KJ * -0.890,
     962             :     KCAL_TO_KJ * -0.890,
     963             :     KCAL_TO_KJ * 2.220,
     964             :     KCAL_TO_KJ * 6.9,
     965             :     1. / (ANG_TO_NM * 3.5),
     966             :     0.199,
     967          10 :   };
     968           5 :   valuemap["CP1"] = {
     969             :     ANG3_TO_NM3 * 25.458,
     970             :     KCAL_TO_KJ * -0.187,
     971             :     KCAL_TO_KJ * -0.187,
     972             :     KCAL_TO_KJ * 0.876,
     973             :     KCAL_TO_KJ * 0.0,
     974             :     1. / (ANG_TO_NM * 3.5),
     975             :     0.227,
     976          10 :   };
     977           5 :   valuemap["CP2"] = {
     978             :     ANG3_TO_NM3 * 19.880,
     979             :     KCAL_TO_KJ * 0.372,
     980             :     KCAL_TO_KJ * 0.372,
     981             :     KCAL_TO_KJ * -0.610,
     982             :     KCAL_TO_KJ * 18.6,
     983             :     1. / (ANG_TO_NM * 3.5),
     984             :     0.217,
     985          10 :   };
     986           5 :   valuemap["CP3"] = {
     987             :     ANG3_TO_NM3 * 26.731,
     988             :     KCAL_TO_KJ * 0.372,
     989             :     KCAL_TO_KJ * 0.372,
     990             :     KCAL_TO_KJ * -0.610,
     991             :     KCAL_TO_KJ * 18.6,
     992             :     1. / (ANG_TO_NM * 3.5),
     993             :     0.217,
     994          10 :   };
     995           5 :   valuemap["CC"] = {
     996             :     ANG3_TO_NM3 * 16.539,
     997             :     KCAL_TO_KJ * 0.000,
     998             :     KCAL_TO_KJ * 0.000,
     999             :     KCAL_TO_KJ * 0.000,
    1000             :     KCAL_TO_KJ * 0.0,
    1001             :     1. / (ANG_TO_NM * 3.5),
    1002             :     0.20,
    1003          10 :   };
    1004           5 :   valuemap["CAI"] = {
    1005             :     ANG3_TO_NM3 * 18.249,
    1006             :     KCAL_TO_KJ * 0.057,
    1007             :     KCAL_TO_KJ * 0.057,
    1008             :     KCAL_TO_KJ * -0.973,
    1009             :     KCAL_TO_KJ * 6.9,
    1010             :     1. / (ANG_TO_NM * 3.5),
    1011             :     0.199,
    1012          10 :   };
    1013           5 :   valuemap["CA"] = {
    1014             :     ANG3_TO_NM3 * 18.249,
    1015             :     KCAL_TO_KJ * 0.057,
    1016             :     KCAL_TO_KJ * 0.057,
    1017             :     KCAL_TO_KJ * -0.973,
    1018             :     KCAL_TO_KJ * 6.9,
    1019             :     1. / (ANG_TO_NM * 3.5),
    1020             :     0.199,
    1021          10 :   };
    1022           5 :   valuemap["N"] = {
    1023             :     ANG3_TO_NM3 * 0.000,
    1024             :     KCAL_TO_KJ * -1.000,
    1025             :     KCAL_TO_KJ * -1.000,
    1026             :     KCAL_TO_KJ * -1.250,
    1027             :     KCAL_TO_KJ * 8.8,
    1028             :     1. / (ANG_TO_NM * 3.5),
    1029             :     0.185,
    1030          10 :   };
    1031           5 :   valuemap["NR1"] = {
    1032             :     ANG3_TO_NM3 * 15.273,
    1033             :     KCAL_TO_KJ * -5.950,
    1034             :     KCAL_TO_KJ * -5.950,
    1035             :     KCAL_TO_KJ * -9.059,
    1036             :     KCAL_TO_KJ * -8.8,
    1037             :     1. / (ANG_TO_NM * 3.5),
    1038             :     0.185,
    1039          10 :   };
    1040           5 :   valuemap["NR2"] = {
    1041             :     ANG3_TO_NM3 * 15.111,
    1042             :     KCAL_TO_KJ * -3.820,
    1043             :     KCAL_TO_KJ * -3.820,
    1044             :     KCAL_TO_KJ * -4.654,
    1045             :     KCAL_TO_KJ * -8.8,
    1046             :     1. / (ANG_TO_NM * 3.5),
    1047             :     0.185,
    1048          10 :   };
    1049           5 :   valuemap["NR3"] = {
    1050             :     ANG3_TO_NM3 * 15.071,
    1051             :     KCAL_TO_KJ * -5.950,
    1052             :     KCAL_TO_KJ * -5.950,
    1053             :     KCAL_TO_KJ * -9.059,
    1054             :     KCAL_TO_KJ * -8.8,
    1055             :     1. / (ANG_TO_NM * 3.5),
    1056             :     0.185,
    1057          10 :   };
    1058           5 :   valuemap["NH1"] = {
    1059             :     ANG3_TO_NM3 * 10.197,
    1060             :     KCAL_TO_KJ * -5.950,
    1061             :     KCAL_TO_KJ * -5.950,
    1062             :     KCAL_TO_KJ * -9.059,
    1063             :     KCAL_TO_KJ * -8.8,
    1064             :     1. / (ANG_TO_NM * 3.5),
    1065             :     0.185,
    1066          10 :   };
    1067           5 :   valuemap["NH2"] = {
    1068             :     ANG3_TO_NM3 * 18.182,
    1069             :     KCAL_TO_KJ * -5.950,
    1070             :     KCAL_TO_KJ * -5.950,
    1071             :     KCAL_TO_KJ * -9.059,
    1072             :     KCAL_TO_KJ * -8.8,
    1073             :     1. / (ANG_TO_NM * 3.5),
    1074             :     0.185,
    1075          10 :   };
    1076           5 :   valuemap["NH3"] = {
    1077             :     ANG3_TO_NM3 * 18.817,
    1078             :     KCAL_TO_KJ * -20.000,
    1079             :     KCAL_TO_KJ * -20.000,
    1080             :     KCAL_TO_KJ * -25.000,
    1081             :     KCAL_TO_KJ * -18.0,
    1082             :     1. / (ANG_TO_NM * 6.0),
    1083             :     0.185,
    1084          10 :   };
    1085           5 :   valuemap["NC2"] = {
    1086             :     ANG3_TO_NM3 * 18.215,
    1087             :     KCAL_TO_KJ * -10.000,
    1088             :     KCAL_TO_KJ * -10.000,
    1089             :     KCAL_TO_KJ * -12.000,
    1090             :     KCAL_TO_KJ * -7.0,
    1091             :     1. / (ANG_TO_NM * 6.0),
    1092             :     0.185,
    1093          10 :   };
    1094           5 :   valuemap["NY"] = {
    1095             :     ANG3_TO_NM3 * 12.001,
    1096             :     KCAL_TO_KJ * -5.950,
    1097             :     KCAL_TO_KJ * -5.950,
    1098             :     KCAL_TO_KJ * -9.059,
    1099             :     KCAL_TO_KJ * -8.8,
    1100             :     1. / (ANG_TO_NM * 3.5),
    1101             :     0.185,
    1102          10 :   };
    1103           5 :   valuemap["NP"] = {
    1104             :     ANG3_TO_NM3 * 4.993,
    1105             :     KCAL_TO_KJ * -20.000,
    1106             :     KCAL_TO_KJ * -20.000,
    1107             :     KCAL_TO_KJ * -25.000,
    1108             :     KCAL_TO_KJ * -18.0,
    1109             :     1. / (ANG_TO_NM * 6.0),
    1110             :     0.185,
    1111          10 :   };
    1112           5 :   valuemap["O"] = {
    1113             :     ANG3_TO_NM3 * 11.772,
    1114             :     KCAL_TO_KJ * -5.330,
    1115             :     KCAL_TO_KJ * -5.330,
    1116             :     KCAL_TO_KJ * -5.787,
    1117             :     KCAL_TO_KJ * -8.8,
    1118             :     1. / (ANG_TO_NM * 3.5),
    1119             :     0.170,
    1120          10 :   };
    1121           5 :   valuemap["OB"] = {
    1122             :     ANG3_TO_NM3 * 11.694,
    1123             :     KCAL_TO_KJ * -5.330,
    1124             :     KCAL_TO_KJ * -5.330,
    1125             :     KCAL_TO_KJ * -5.787,
    1126             :     KCAL_TO_KJ * -8.8,
    1127             :     1. / (ANG_TO_NM * 3.5),
    1128             :     0.170,
    1129          10 :   };
    1130           5 :   valuemap["OC"] = {
    1131             :     ANG3_TO_NM3 * 12.003,
    1132             :     KCAL_TO_KJ * -10.000,
    1133             :     KCAL_TO_KJ * -10.000,
    1134             :     KCAL_TO_KJ * -12.000,
    1135             :     KCAL_TO_KJ * -9.4,
    1136             :     1. / (ANG_TO_NM * 6.0),
    1137             :     0.170,
    1138          10 :   };
    1139           5 :   valuemap["OH1"] = {
    1140             :     ANG3_TO_NM3 * 15.528,
    1141             :     KCAL_TO_KJ * -5.920,
    1142             :     KCAL_TO_KJ * -5.920,
    1143             :     KCAL_TO_KJ * -9.264,
    1144             :     KCAL_TO_KJ * -11.2,
    1145             :     1. / (ANG_TO_NM * 3.5),
    1146             :     0.177,
    1147          10 :   };
    1148           5 :   valuemap["OS"] = {
    1149             :     ANG3_TO_NM3 * 6.774,
    1150             :     KCAL_TO_KJ * -2.900,
    1151             :     KCAL_TO_KJ * -2.900,
    1152             :     KCAL_TO_KJ * -3.150,
    1153             :     KCAL_TO_KJ * -4.8,
    1154             :     1. / (ANG_TO_NM * 3.5),
    1155             :     0.177,
    1156          10 :   };
    1157           5 :   valuemap["S"] = {
    1158             :     ANG3_TO_NM3 * 20.703,
    1159             :     KCAL_TO_KJ * -3.240,
    1160             :     KCAL_TO_KJ * -3.240,
    1161             :     KCAL_TO_KJ * -4.475,
    1162             :     KCAL_TO_KJ * -39.9,
    1163             :     1. / (ANG_TO_NM * 3.5),
    1164             :     0.20,
    1165          10 :   };
    1166           5 :   valuemap["SM"] = {
    1167             :     ANG3_TO_NM3 * 21.306,
    1168             :     KCAL_TO_KJ * -3.240,
    1169             :     KCAL_TO_KJ * -3.240,
    1170             :     KCAL_TO_KJ * -4.475,
    1171             :     KCAL_TO_KJ * -39.9,
    1172             :     1. / (ANG_TO_NM * 3.5),
    1173             :     0.197,
    1174          10 :   };
    1175           5 :   return valuemap;
    1176             : }
    1177             : }
    1178             : }

Generated by: LCOV version 1.16