Remove ResidueTypeMap dependency from gmx chi
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 7 Sep 2021 05:21:36 +0000 (07:21 +0200)
committerAndrey Alekseenko <al42and@gmail.com>
Wed, 22 Sep 2021 08:35:17 +0000 (08:35 +0000)
This dependency had the sole effect of making the tool harder to use
when residue names unknown to GROMACS were used. The dihedrals are
found from the atom names, and the residue name is used only for
reporting results, so there is no reason to require the user to add
the name to an arbitrary list.

docs/release-notes/2022/major/tools.rst
src/gromacs/gmxana/dlist.cpp
src/gromacs/gmxana/gmx_chi.cpp
src/gromacs/gmxana/gstat.h

index fdc5086b6992a71ceacef0503ddb83590d5ae4de..e35b1f3e8f05bcb4355399eb54dfe051d64d365d 100644 (file)
@@ -7,8 +7,8 @@ Improvements to |Gromacs| tools
    Also, please use the syntax :issue:`number` to reference issues on GitLab, without the
    a space between the colon and number!
 
-gmx msd has been migrated to the trajectoryanalysis framework
-"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+``gmx msd`` has been migrated to the trajectoryanalysis framework
+"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
 
 The tool now uses the |Gromacs| selection syntax. Rather than piping selections via stdin,
 selections are now made using the "-sel" option.
@@ -23,4 +23,10 @@ Some rarely used features have yet to be migrated, including:
 - System COM removal with -rmcomm has not yet been implemented.
 - B-factor writing using the -pdb option is not yet supported.
 
-:issue:`2368`
\ No newline at end of file
+:issue:`2368`
+
+``gmx chi`` no longer needs ``residuetypes.dat`` entries for custom residues
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+The need to add the names of custom residues to ``residuetypes.dat`` has been
+removed, because it served no purpose. This makes ``gmx chi`` easier to use.
index c6536b8aa56fcd2bb25e217c429336563ab0c464..d4fe29fc644e02971f6d1abdc1ac0abafc0875e6 100644 (file)
 #include "gromacs/topology/topology.h"
 #include "gromacs/utility/fatalerror.h"
 
-std::vector<t_dlist> mk_dlist(FILE*                 log,
-                              const t_atoms*        atoms,
-                              gmx_bool              bPhi,
-                              gmx_bool              bPsi,
-                              gmx_bool              bChi,
-                              gmx_bool              bHChi,
-                              int                   maxchi,
-                              int                   r0,
-                              const ResidueTypeMap& residueTypeMap)
+std::vector<t_dlist> mk_dlist(FILE*          log,
+                              const t_atoms* atoms,
+                              gmx_bool       bPhi,
+                              gmx_bool       bPsi,
+                              gmx_bool       bChi,
+                              gmx_bool       bHChi,
+                              int            maxchi,
+                              int            r0)
 {
     int       i, j, ii;
     t_dihatms atm, prev;
@@ -227,16 +226,6 @@ std::vector<t_dlist> mk_dlist(FILE*                 log,
                 nc[6]++;
             }
 
-            /* Prevent use of unknown residues. If one adds a custom residue to
-             * residuetypes.dat but somehow loses it, changes it, or does analysis on
-             * another machine, the residue type will be unknown. */
-            if (residueTypeMap.find(thisres) == residueTypeMap.end())
-            {
-                gmx_fatal(FARGS,
-                          "Unknown residue %s when searching for residue type.\n"
-                          "Maybe you need to add a custom residue in residuetypes.dat.",
-                          thisres);
-            }
             dl[nl].residueName = thisres;
 
             sprintf(dl[nl].name, "%s%d", thisres, ires + r0);
index 0ca1075c6e02d7be8446d0c10c10a58a6d130cc9..515224cd04f857b2e53bd3c5256b20f7812d0f21 100644 (file)
@@ -60,7 +60,6 @@
 #include "gromacs/math/units.h"
 #include "gromacs/math/utilities.h"
 #include "gromacs/math/vec.h"
-#include "gromacs/topology/residuetypes.h"
 #include "gromacs/topology/topology.h"
 #include "gromacs/utility/arraysize.h"
 #include "gromacs/utility/cstringutil.h"
@@ -519,10 +518,7 @@ static void histogramming(FILE*                    log,
     }
 
     // Build a list of unique residue names found in the dihedral
-    // list, so we can loop over those unique names conveniently. The
-    // names are the same as the residue names found in residueTypeMap in the
-    // caller, but ResidueTypeMap doesn't yet have a way to loop over its
-    // contents.
+    // list, so we can loop over those unique names conveniently.
     std::unordered_set<std::string> uniqueResidueNames;
     for (const auto& dihedral : dlist)
     {
@@ -1530,9 +1526,7 @@ int gmx_chi(int argc, char* argv[])
     }
     fprintf(log, "Title: %s\n", name);
 
-    ResidueTypeMap       residueTypeMap = residueTypeMapFromLibraryFile("residuetypes.dat");
-    std::vector<t_dlist> dlist =
-            mk_dlist(log, &atoms, bPhi, bPsi, bChi, bHChi, maxchi, r0, residueTypeMap);
+    std::vector<t_dlist> dlist = mk_dlist(log, &atoms, bPhi, bPsi, bChi, bHChi, maxchi, r0);
     fprintf(stderr, "%zu residues with dihedrals found\n", dlist.size());
 
     if (dlist.empty())
index 32075660630b1d330457df19a36d89dec5e9711c..ec13293808ef17b45e4760fa881e027d8be6e97d 100644 (file)
@@ -42,7 +42,6 @@
 
 #include "gromacs/commandline/pargs.h"
 #include "gromacs/topology/index.h"
-#include "gromacs/topology/residuetypes.h"
 
 struct gmx_output_env_t;
 
@@ -255,15 +254,14 @@ gmx_bool has_dihedral(int Dih, const t_dlist& dlist);
  * the residues, and a mapping from chemical peptide atom names to
  * atom indices based on the atom names. Many fields of t_dlist are
  * not yet filled. */
-std::vector<t_dlist> mk_dlist(FILE*                 log,
-                              const t_atoms*        atoms,
-                              gmx_bool              bPhi,
-                              gmx_bool              bPsi,
-                              gmx_bool              bChi,
-                              gmx_bool              bHChi,
-                              int                   maxchi,
-                              int                   r0,
-                              const ResidueTypeMap& rt);
+std::vector<t_dlist> mk_dlist(FILE*          log,
+                              const t_atoms* atoms,
+                              gmx_bool       bPhi,
+                              gmx_bool       bPsi,
+                              gmx_bool       bChi,
+                              gmx_bool       bHChi,
+                              int            maxchi,
+                              int            r0);
 
 void pr_dlist(FILE*                        fp,
               gmx::ArrayRef<const t_dlist> dlist,