LCOV - code coverage report
Current view: top level - colvar - EEFSolv.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 149 166 89.8 %
Date: 2021-11-18 15:22:58 Functions: 15 16 93.8 %

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

Generated by: LCOV version 1.14