Introduce gmx report-methods
authorPaul Bauer <paul.bauer.q@gmail.com>
Thu, 17 May 2018 18:04:57 +0000 (20:04 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Fri, 14 Sep 2018 08:54:45 +0000 (10:54 +0200)
The writing of system information from gmx check has been
moved to a new module gmx report-methods that does the only
this part, writing the basic system information to the
the terminal, as well as a LaTeX formatted or unformatted
output file.

Change-Id: Id34c8ca67e606a88da70fe9558b2811a716dcc47

src/gromacs/tools/CMakeLists.txt
src/gromacs/tools/report-methods.cpp [new file with mode: 0644]
src/gromacs/tools/report-methods.h [new file with mode: 0644]
src/gromacs/tools/tests/CMakeLists.txt [new file with mode: 0644]
src/gromacs/tools/tests/lysozyme.pdb [new file with mode: 0644]
src/gromacs/tools/tests/lysozyme.top [new file with mode: 0644]
src/gromacs/tools/tests/report-methods.cpp [new file with mode: 0644]
src/programs/legacymodules.cpp

index 061fcf2168067e3762951214c67d1e13c90dc179..442b77b234837fa9b31f0f6d6cac5abb958577e6 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2014,2015, by the GROMACS development team, led by
+# Copyright (c) 2014,2015,2018, by the GROMACS development team, led by
 # Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
 # and including many others, as listed in the AUTHORS file in the
 # top-level source directory and at http://www.gromacs.org.
@@ -36,5 +36,5 @@ file(GLOB TOOLS_SOURCES *.cpp)
 set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${TOOLS_SOURCES} PARENT_SCOPE)
 
 if (BUILD_TESTING)
-#    add_subdirectory(tests)
+    add_subdirectory(tests)
 endif()
diff --git a/src/gromacs/tools/report-methods.cpp b/src/gromacs/tools/report-methods.cpp
new file mode 100644 (file)
index 0000000..3700c3d
--- /dev/null
@@ -0,0 +1,247 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+#include "gmxpre.h"
+
+#include "report-methods.h"
+
+#include "gromacs/commandline/cmdlineoptionsmodule.h"
+#include "gromacs/fileio/confio.h"
+#include "gromacs/fileio/filetypes.h"
+#include "gromacs/fileio/tpxio.h"
+#include "gromacs/mdtypes/inputrec.h"
+#include "gromacs/mdtypes/state.h"
+#include "gromacs/options/basicoptions.h"
+#include "gromacs/options/filenameoption.h"
+#include "gromacs/options/ioptionscontainer.h"
+#include "gromacs/selection/selectionoptionbehavior.h"
+#include "gromacs/topology/mtop_util.h"
+#include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/fileredirector.h"
+#include "gromacs/utility/filestream.h"
+#include "gromacs/utility/stringutil.h"
+
+namespace gmx
+{
+
+void writeHeader(TextWriter *writer, const std::string &text, const std::string &section, bool writeFormattedText)
+{
+    std::string formattedText;
+    if (writeFormattedText)
+    {
+        formattedText = "\\" + section + "{" + text + "}\n";
+    }
+    else
+    {
+        formattedText = section + ": " + text + "\n";
+    }
+    writer->writeString(formattedText);
+}
+
+void writeSystemInformation(TextWriter *writer, const gmx_mtop_t &top, bool writeFormattedText)
+{
+    int                       nmol, nvsite = 0;
+    gmx_mtop_atomloop_block_t aloop;
+    const t_atom             *atom;
+
+    writeHeader(writer, "Simulation system", "subsection", writeFormattedText);
+    aloop = gmx_mtop_atomloop_block_init(&top);
+    while (gmx_mtop_atomloop_block_next(aloop, &atom, &nmol))
+    {
+        if (atom->ptype == eptVSite)
+        {
+            nvsite += nmol;
+        }
+    }
+    {
+        writer->writeLine(formatString("A system of %d molecules (%d atoms) was simulated.",
+                                       gmx_mtop_num_molecules(top), top.natoms-nvsite));
+    }
+    if (nvsite)
+    {
+        writer->writeLine(formatString("Virtual sites were used in some of the molecules."));
+    }
+    writer->ensureEmptyLine();
+}
+
+void writeParameterInformation(TextWriter *writer, const t_inputrec &ir, bool writeFormattedText)
+{
+    writeHeader(writer, "Simulation settings", "subsection", writeFormattedText);
+    writer->writeLine(formatString("A total of %g ns were simulated with a time step of %g fs.",
+                                   ir.nsteps*ir.delta_t*0.001, 1000*ir.delta_t));
+    writer->writeLine(formatString("Neighbor searching was performed every %d steps.", ir.nstlist));
+    writer->writeLine(formatString("The %s algorithm was used for electrostatic interactions.",
+                                   EELTYPE(ir.coulombtype)));
+    writer->writeLine(formatString("with a cut-off of %g nm.", ir.rcoulomb));
+    if (ir.coulombtype == eelPME)
+    {
+        writer->writeLine(formatString("A reciprocal grid of %d x %d x %d cells was used with %dth order B-spline interpolation.", ir.nkx, ir.nky, ir.nkz, ir.pme_order));
+    }
+    writer->writeLine(formatString("A single cut-off of %g nm was used for Van der Waals interactions.", ir.rlist));
+    if (ir.etc != 0)
+    {
+        writer->writeLine(formatString("Temperature coupling was done with the %s algorithm.",
+                                       etcoupl_names[ir.etc]));
+    }
+    if (ir.epc != 0)
+    {
+        writer->writeLine(formatString("Pressure coupling was done with the %s algorithm.",
+                                       epcoupl_names[ir.epc]));
+    }
+    writer->ensureEmptyLine();
+}
+
+void writeInformation(TextOutputFile *outputStream, const t_inputrec &ir,
+                      const gmx_mtop_t &top, bool writeFormattedText, bool notStdout)
+{
+    TextWriter              writer(outputStream);
+    writer.ensureEmptyLine();
+    writeHeader(&writer, "Methods", "section", writeFormattedText);
+    writeSystemInformation(&writer, top, writeFormattedText);
+    writeParameterInformation(&writer, ir, writeFormattedText);
+    writer.ensureEmptyLine();
+
+    if (notStdout)
+    {
+        writer.close();
+    }
+}
+
+namespace
+{
+
+class ReportMethods : public ICommandLineOptionsModule
+{
+    public:
+        ReportMethods()
+            : writeLatex_(false), writePlainText_(false)
+        {
+        }
+
+        // From ICommandLineOptionsModule
+        void init(CommandLineModuleSettings * /*settings*/) override
+        {
+        }
+        void initOptions(IOptionsContainer                 *options,
+                         ICommandLineOptionsModuleSettings *settings) override;
+        void optionsFinished() override;
+        int run() override;
+
+    private:
+
+        //! File name for the output LaTeX file or empty.
+        std::string         outputFileLatex_;
+        //! File name for the unformatted output file or empty.
+        std::string         outputFileUnformatted_;
+        //! File name of the run input file with full topology.
+        std::string         inputTopology_;
+        //! Boolean reporting if writing to the LaTeX output file is requested.
+        bool                writeLatex_;
+        //! Boolean reporting if writing to unformatted output is requested.
+        bool                writePlainText_;
+};
+
+void ReportMethods::initOptions(IOptionsContainer                 *options,
+                                ICommandLineOptionsModuleSettings *settings)
+{
+    const char *const desc[] = {
+        "[THISMODULE] reports basic system information for the run input",
+        "file specfied with [TT]-s[tt] either to the",
+        "terminal, to a LaTeX formatted output file if run with",
+        "the [TT]-m[tt] option or to an unformatted file with",
+        "the [TT]-o[tt] option.",
+        "The functionality has been moved here from its previous",
+        "place in [gmx-check]."
+    };
+
+    settings->setHelpText(desc);
+
+    options->addOption(FileNameOption("s")
+                           .filetype(eftTopology).inputFile().required()
+                           .store(&inputTopology_)
+                           .defaultBasename("topol")
+                           .description("Run input file for report"));
+
+    // TODO: Replace use of legacyType.
+    options->addOption(FileNameOption("m")
+                           .legacyType(efTEX).outputFile()
+                           .store(&outputFileLatex_).storeIsSet(&writeLatex_)
+                           .defaultBasename("report")
+                           .description("LaTeX formatted report output"));
+    options->addOption(FileNameOption("o")
+                           .legacyType(efOUT).outputFile()
+                           .store(&outputFileUnformatted_).storeIsSet(&writePlainText_)
+                           .defaultBasename("report")
+                           .description("Unformatted report output to file"));
+
+}
+
+void ReportMethods::optionsFinished()
+{
+}
+
+int ReportMethods::run()
+{
+    t_state    state;
+    t_inputrec ir;
+    gmx_mtop_t top;
+    read_tpx_state(inputTopology_.c_str(), &ir, &state, &top);
+    if (writeLatex_)
+    {
+        TextOutputFile file(outputFileLatex_);
+        writeInformation(&file, ir, top, true, true);
+    }
+    if (writePlainText_)
+    {
+        TextOutputFile file(outputFileUnformatted_);
+        writeInformation(&file, ir, top, false, true);
+    }
+    TextOutputFile         &stdoutFile = TextOutputFile::standardOutput();
+    writeInformation(&stdoutFile, ir, top, false, false);
+
+    return 0;
+}
+
+}   // namespace
+
+const char ReportMethodsInfo::name[]             = "report";
+const char ReportMethodsInfo::shortDescription[] =
+    "Write short summary about the simulation setup to a text file "
+    "and/or to the standard output.";
+ICommandLineOptionsModulePointer ReportMethodsInfo::create()
+{
+    return ICommandLineOptionsModulePointer(compat::make_unique<ReportMethods>());
+}
+
+} // namespace gmx
diff --git a/src/gromacs/tools/report-methods.h b/src/gromacs/tools/report-methods.h
new file mode 100644 (file)
index 0000000..221cdf2
--- /dev/null
@@ -0,0 +1,110 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+#ifndef GMX_TOOLS_REPORT_METHODS_H
+#define GMX_TOOLS_REPORT_METHODS_H
+
+#include "gromacs/commandline/cmdlineoptionsmodule.h"
+#include "gromacs/mdtypes/inputrec.h"
+#include "gromacs/topology/topology.h"
+#include "gromacs/utility/filestream.h"
+#include "gromacs/utility/textwriter.h"
+
+namespace gmx
+{
+
+class ReportMethodsInfo
+{
+    public:
+        static const char name[];
+        static const char shortDescription[];
+        static ICommandLineOptionsModulePointer create();
+};
+
+// Helper functions of the class
+
+/*! \brief
+ * Write appropiate Header to output stream.
+ *
+ * \param[in] writer TextWriter object for writing information.
+ * \param[in] text String with the header before writing.
+ * \param[in] section String with section text for header.
+ * \param[in] writeFormattedText If we need to format the text for LaTeX output or not
+ */
+void writeHeader(TextWriter *writer, const std::string &text, const std::string &section, bool writeFormattedText);
+
+/*! \brief
+ * Write information about the molecules in the system.
+ *
+ * This method should write all possible information about
+ * the molecular composition of the system.
+ *
+ * \param[in] writer TextWriter object for writing information.
+ * \param[in] top Local topology used to derive the information to write out.
+ * \param[in] writeFormattedText Decide if we want formatted text output or not.
+ *
+ */
+void writeSystemInformation(TextWriter *writer, const gmx_mtop_t &top, bool writeFormattedText);
+
+/*! \brief
+ * Write information about system parameters.
+ *
+ * This method writes the basic information for the system parameters
+ * and simulation settings as reported in the \p ir.
+ *
+ * \param[in] writer TextWriter object for writing information.
+ * \param[in] ir Reference to inputrec of the run input.
+ * \param[in] writeFormattedText Decide if we want formatted text output or not.
+ */
+void writeParameterInformation(TextWriter *writer, const t_inputrec &ir, bool writeFormattedText);
+
+/*! \brief
+ * Wrapper for writing out information.
+ *
+ * This function is actual called from within the run method
+ * to write the information to the terminal or to file.
+ * New write out methods should be added to it instead of adding them in run.
+ *
+ * \param[in] outputStream The filestream used to write the information to.
+ * \param[in] ir Reference to inputrec of the run input.
+ * \param[in] top Local topology used to derive the information to write out.
+ * \param[in] writeFormattedText Decide if we want formatted text output or not.
+ * \param[in] notStdout Bool to see if we can close the file after writing or not in case of stdout.
+ */
+void writeInformation(TextOutputFile *outputStream, const t_inputrec &ir,
+                      const gmx_mtop_t &top, bool writeFormattedText, bool notStdout);
+
+} // namespace gmx
+
+#endif
diff --git a/src/gromacs/tools/tests/CMakeLists.txt b/src/gromacs/tools/tests/CMakeLists.txt
new file mode 100644 (file)
index 0000000..d543a8c
--- /dev/null
@@ -0,0 +1,37 @@
+#
+# This file is part of the GROMACS molecular simulation package.
+#
+# Copyright (c) 2018, by the GROMACS development team, led by
+# Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+# and including many others, as listed in the AUTHORS file in the
+# top-level source directory and at http://www.gromacs.org.
+#
+# GROMACS is free software; you can redistribute it and/or
+# modify it under the terms of the GNU Lesser General Public License
+# as published by the Free Software Foundation; either version 2.1
+# of the License, or (at your option) any later version.
+#
+# GROMACS is distributed in the hope that it will be useful,
+# but WITHOUT ANY WARRANTY; without even the implied warranty of
+# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+# Lesser General Public License for more details.
+#
+# You should have received a copy of the GNU Lesser General Public
+# License along with GROMACS; if not, see
+# http://www.gnu.org/licenses, or write to the Free Software Foundation,
+# Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+#
+# If you want to redistribute modifications to GROMACS, please
+# consider that scientific software is very special. Version
+# control is crucial - bugs must be traceable. We will be happy to
+# consider code for inclusion in the official distribution, but
+# derived work must not be called official GROMACS. Details are found
+# in the README & COPYING files - if they are missing, get the
+# official version at http://www.gromacs.org.
+#
+# To help us fund GROMACS development, we humbly ask that you cite
+# the research papers on the package. Check out http://www.gromacs.org.
+
+gmx_add_unit_test(ToolUnitTests tool-test
+                  report-methods.cpp)
+
diff --git a/src/gromacs/tools/tests/lysozyme.pdb b/src/gromacs/tools/tests/lysozyme.pdb
new file mode 100644 (file)
index 0000000..d82163b
--- /dev/null
@@ -0,0 +1,162 @@
+TITLE     First 10 residues from 1AKI
+REMARK    THIS IS A SIMULATION BOX
+CRYST1   59.062   68.451   30.517  90.00  90.00  90.00 P 1           1
+MODEL        1
+ATOM      1  N   LYS B   1      35.360  22.340 -11.980  1.00  0.00           N
+ATOM      2  H1  LYS B   1      36.120  22.880 -12.360  1.00  0.00           H
+ATOM      3  H2  LYS B   1      34.700  22.140 -12.700  1.00  0.00           H
+ATOM      4  H3  LYS B   1      34.920  22.860 -11.250  1.00  0.00           H
+ATOM      5  CA  LYS B   1      35.890  21.070 -11.430  1.00  0.00           C
+ATOM      6  HA  LYS B   1      36.330  20.550 -12.160  1.00  0.00           H
+ATOM      7  CB  LYS B   1      36.870  21.440 -10.310  1.00  0.00           C
+ATOM      8  HB1 LYS B   1      37.630  21.950 -10.700  1.00  0.00           H
+ATOM      9  HB2 LYS B   1      36.390  22.010  -9.640  1.00  0.00           H
+ATOM     10  CG  LYS B   1      37.450  20.250  -9.560  1.00  0.00           C
+ATOM     11  HG1 LYS B   1      36.760  19.890  -8.940  1.00  0.00           H
+ATOM     12  HG2 LYS B   1      37.700  19.540 -10.230  1.00  0.00           H
+ATOM     13  CD  LYS B   1      38.690  20.650  -8.770  1.00  0.00           C
+ATOM     14  HD1 LYS B   1      39.450  20.830  -9.400  1.00  0.00           H
+ATOM     15  HD2 LYS B   1      38.490  21.470  -8.240  1.00  0.00           H
+ATOM     16  CE  LYS B   1      39.060  19.510  -7.840  1.00  0.00           C
+ATOM     17  HE1 LYS B   1      38.410  19.460  -7.080  1.00  0.00           H
+ATOM     18  HE2 LYS B   1      39.060  18.640  -8.330  1.00  0.00           H
+ATOM     19  NZ  LYS B   1      40.420  19.770  -7.300  1.00  0.00           N
+ATOM     20  HZ1 LYS B   1      40.690  19.030  -6.680  1.00  0.00           H
+ATOM     21  HZ2 LYS B   1      41.080  19.820  -8.060  1.00  0.00           H
+ATOM     22  HZ3 LYS B   1      40.420  20.640  -6.800  1.00  0.00           H
+ATOM     23  C   LYS B   1      34.740  20.260 -10.840  1.00  0.00           C
+ATOM     24  O   LYS B   1      33.950  20.810 -10.080  1.00  0.00           O
+ATOM     25  N   VAL B   2      34.740  18.960 -11.040  1.00  0.00           N
+ATOM     26  H   VAL B   2      35.360  18.600 -11.740  1.00  0.00           H
+ATOM     27  CA  VAL B   2      33.900  18.000 -10.330  1.00  0.00           C
+ATOM     28  HA  VAL B   2      33.170  18.520  -9.900  1.00  0.00           H
+ATOM     29  CB  VAL B   2      33.140  17.030 -11.230  1.00  0.00           C
+ATOM     30  HB  VAL B   2      33.860  16.520 -11.700  1.00  0.00           H
+ATOM     31  CG1 VAL B   2      32.250  16.080 -10.430  1.00  0.00           C
+ATOM     32 1HG1 VAL B   2      31.770  15.470 -11.060  1.00  0.00           H
+ATOM     33 2HG1 VAL B   2      32.820  15.550  -9.810  1.00  0.00           H
+ATOM     34 3HG1 VAL B   2      31.580  16.610  -9.910  1.00  0.00           H
+ATOM     35  CG2 VAL B   2      32.290  17.710 -12.290  1.00  0.00           C
+ATOM     36 1HG2 VAL B   2      31.830  17.020 -12.840  1.00  0.00           H
+ATOM     37 2HG2 VAL B   2      31.620  18.300 -11.850  1.00  0.00           H
+ATOM     38 3HG2 VAL B   2      32.880  18.270 -12.880  1.00  0.00           H
+ATOM     39  C   VAL B   2      34.800  17.310  -9.290  1.00  0.00           C
+ATOM     40  O   VAL B   2      35.760  16.610  -9.660  1.00  0.00           O
+ATOM     41  N   PHE B   3      34.490  17.550  -8.040  1.00  0.00           N
+ATOM     42  H   PHE B   3      33.750  18.190  -7.840  1.00  0.00           H
+ATOM     43  CA  PHE B   3      35.190  16.900  -6.920  1.00  0.00           C
+ATOM     44  HA  PHE B   3      36.150  16.970  -7.170  1.00  0.00           H
+ATOM     45  CB  PHE B   3      34.970  17.630  -5.590  1.00  0.00           C
+ATOM     46  HB1 PHE B   3      34.050  18.020  -5.580  1.00  0.00           H
+ATOM     47  HB2 PHE B   3      35.060  16.980  -4.840  1.00  0.00           H
+ATOM     48  CG  PHE B   3      35.940  18.740  -5.380  1.00  0.00           C
+ATOM     49  CD1 PHE B   3      35.670  20.050  -5.800  1.00  0.00           C
+ATOM     50  HD1 PHE B   3      34.810  20.250  -6.270  1.00  0.00           H
+ATOM     51  CD2 PHE B   3      37.000  18.560  -4.470  1.00  0.00           C
+ATOM     52  HD2 PHE B   3      37.130  17.660  -4.050  1.00  0.00           H
+ATOM     53  CE1 PHE B   3      36.580  21.080  -5.570  1.00  0.00           C
+ATOM     54  HE1 PHE B   3      36.480  21.950  -6.040  1.00  0.00           H
+ATOM     55  CE2 PHE B   3      37.870  19.590  -4.160  1.00  0.00           C
+ATOM     56  HE2 PHE B   3      38.660  19.420  -3.570  1.00  0.00           H
+ATOM     57  CZ  PHE B   3      37.640  20.870  -4.670  1.00  0.00           C
+ATOM     58  HZ  PHE B   3      38.220  21.640  -4.390  1.00  0.00           H
+ATOM     59  C   PHE B   3      34.740  15.440  -6.770  1.00  0.00           C
+ATOM     60  O   PHE B   3      33.520  15.160  -6.860  1.00  0.00           O
+ATOM     61  N   GLY B   4      35.720  14.640  -6.330  1.00  0.00           N
+ATOM     62  H   GLY B   4      36.670  14.950  -6.320  1.00  0.00           H
+ATOM     63  CA  GLY B   4      35.370  13.280  -5.870  1.00  0.00           C
+ATOM     64  HA1 GLY B   4      34.620  12.920  -6.430  1.00  0.00           H
+ATOM     65  HA2 GLY B   4      36.160  12.680  -5.940  1.00  0.00           H
+ATOM     66  C   GLY B   4      34.920  13.420  -4.420  1.00  0.00           C
+ATOM     67  O   GLY B   4      35.300  14.400  -3.780  1.00  0.00           O
+ATOM     68  N   ARG B   5      34.050  12.540  -3.970  1.00  0.00           N
+ATOM     69  H   ARG B   5      33.710  11.840  -4.600  1.00  0.00           H
+ATOM     70  CA  ARG B   5      33.560  12.540  -2.590  1.00  0.00           C
+ATOM     71  HA  ARG B   5      32.980  13.340  -2.520  1.00  0.00           H
+ATOM     72  CB  ARG B   5      32.760  11.260  -2.330  1.00  0.00           C
+ATOM     73  HB1 ARG B   5      32.000  11.220  -2.970  1.00  0.00           H
+ATOM     74  HB2 ARG B   5      33.360  10.470  -2.470  1.00  0.00           H
+ATOM     75  CG  ARG B   5      32.210  11.200  -0.920  1.00  0.00           C
+ATOM     76  HG1 ARG B   5      32.970  11.170  -0.270  1.00  0.00           H
+ATOM     77  HG2 ARG B   5      31.650  12.010  -0.750  1.00  0.00           H
+ATOM     78  CD  ARG B   5      31.380  10.000  -0.720  1.00  0.00           C
+ATOM     79  HD1 ARG B   5      31.040   9.990   0.220  1.00  0.00           H
+ATOM     80  HD2 ARG B   5      30.600  10.050  -1.350  1.00  0.00           H
+ATOM     81  NE  ARG B   5      32.060   8.750  -0.960  1.00  0.00           N
+ATOM     82  HE  ARG B   5      32.020   8.400  -1.890  1.00  0.00           H
+ATOM     83  CZ  ARG B   5      32.730   8.010  -0.100  1.00  0.00           C
+ATOM     84  NH1 ARG B   5      32.840   8.330   1.190  1.00  0.00           N
+ATOM     85 1HH1 ARG B   5      32.390   9.160   1.530  1.00  0.00           H
+ATOM     86 2HH1 ARG B   5      33.360   7.750   1.810  1.00  0.00           H
+ATOM     87  NH2 ARG B   5      33.250   6.840  -0.530  1.00  0.00           N
+ATOM     88 1HH2 ARG B   5      33.110   6.550  -1.470  1.00  0.00           H
+ATOM     89 2HH2 ARG B   5      33.760   6.260   0.100  1.00  0.00           H
+ATOM     90  C   ARG B   5      34.670  12.730  -1.560  1.00  0.00           C
+ATOM     91  O   ARG B   5      34.670  13.650  -0.700  1.00  0.00           O
+ATOM     92  N   CYS B   6      35.670  11.850  -1.610  1.00  0.00           N
+ATOM     93  H   CYS B   6      35.670  11.160  -2.330  1.00  0.00           H
+ATOM     94  CA  CYS B   6      36.780  11.870  -0.650  1.00  0.00           C
+ATOM     95  HA  CYS B   6      36.310  12.020   0.220  1.00  0.00           H
+ATOM     96  CB  CYS B   6      37.490  10.530  -0.620  1.00  0.00           C
+ATOM     97  HB1 CYS B   6      37.700  10.340  -1.580  1.00  0.00           H
+ATOM     98  HB2 CYS B   6      38.340  10.720  -0.130  1.00  0.00           H
+ATOM     99  SG  CYS B   6      36.540   9.200   0.140  1.00  0.00           S
+ATOM    100  HG  CYS B   6      37.075   8.355   0.120  1.00  0.00           H
+ATOM    101  C   CYS B   6      37.750  13.050  -0.780  1.00  0.00           C
+ATOM    102  O   CYS B   6      38.150  13.610   0.260  1.00  0.00           O
+ATOM    103  N   GLU B   7      37.860  13.480  -2.020  1.00  0.00           N
+ATOM    104  H   GLU B   7      37.400  13.000  -2.760  1.00  0.00           H
+ATOM    105  CA  GLU B   7      38.680  14.690  -2.310  1.00  0.00           C
+ATOM    106  HA  GLU B   7      39.600  14.550  -1.930  1.00  0.00           H
+ATOM    107  CB  GLU B   7      38.780  14.850  -3.820  1.00  0.00           C
+ATOM    108  HB1 GLU B   7      39.230  14.020  -4.170  1.00  0.00           H
+ATOM    109  HB2 GLU B   7      37.850  14.890  -4.170  1.00  0.00           H
+ATOM    110  CG  GLU B   7      39.540  16.050  -4.380  1.00  0.00           C
+ATOM    111  HG1 GLU B   7      39.130  16.870  -3.990  1.00  0.00           H
+ATOM    112  HG2 GLU B   7      40.490  15.980  -4.070  1.00  0.00           H
+ATOM    113  CD  GLU B   7      39.580  16.240  -5.870  1.00  0.00           C
+ATOM    114  OE1 GLU B   7      38.670  15.640  -6.490  1.00  0.00           O
+ATOM    115  OE2 GLU B   7      40.420  16.950  -6.380  1.00  0.00           O
+ATOM    116  C   GLU B   7      38.050  15.930  -1.660  1.00  0.00           C
+ATOM    117  O   GLU B   7      38.740  16.730  -1.010  1.00  0.00           O
+ATOM    118  N   LEU B   8      36.740  16.050  -1.820  1.00  0.00           N
+ATOM    119  H   LEU B   8      36.260  15.350  -2.350  1.00  0.00           H
+ATOM    120  CA  LEU B   8      35.960  17.160  -1.250  1.00  0.00           C
+ATOM    121  HA  LEU B   8      36.400  18.010  -1.560  1.00  0.00           H
+ATOM    122  CB  LEU B   8      34.530  17.170  -1.810  1.00  0.00           C
+ATOM    123  HB1 LEU B   8      34.570  17.220  -2.810  1.00  0.00           H
+ATOM    124  HB2 LEU B   8      34.060  16.330  -1.530  1.00  0.00           H
+ATOM    125  CG  LEU B   8      33.720  18.350  -1.310  1.00  0.00           C
+ATOM    126  HG  LEU B   8      33.780  18.420  -0.310  1.00  0.00           H
+ATOM    127  CD1 LEU B   8      34.300  19.660  -1.840  1.00  0.00           C
+ATOM    128 1HD1 LEU B   8      33.760  20.430  -1.500  1.00  0.00           H
+ATOM    129 2HD1 LEU B   8      35.240  19.750  -1.530  1.00  0.00           H
+ATOM    130 3HD1 LEU B   8      34.270  19.650  -2.840  1.00  0.00           H
+ATOM    131  CD2 LEU B   8      32.250  18.140  -1.600  1.00  0.00           C
+ATOM    132 1HD2 LEU B   8      31.720  18.930  -1.260  1.00  0.00           H
+ATOM    133 2HD2 LEU B   8      32.110  18.050  -2.580  1.00  0.00           H
+ATOM    134 3HD2 LEU B   8      31.930  17.310  -1.140  1.00  0.00           H
+ATOM    135  C   LEU B   8      36.050  17.130   0.270  1.00  0.00           C
+ATOM    136  O   LEU B   8      36.160  18.170   0.920  1.00  0.00           O
+ATOM    137  N   ALA B   9      35.750  15.980   0.830  1.00  0.00           N
+ATOM    138  H   ALA B   9      35.460  15.220   0.240  1.00  0.00           H
+ATOM    139  CA  ALA B   9      35.840  15.760   2.280  1.00  0.00           C
+ATOM    140  HA  ALA B   9      35.080  16.260   2.690  1.00  0.00           H
+ATOM    141  CB  ALA B   9      35.660  14.290   2.620  1.00  0.00           C
+ATOM    142  HB1 ALA B   9      35.720  14.160   3.610  1.00  0.00           H
+ATOM    143  HB2 ALA B   9      34.760  13.980   2.300  1.00  0.00           H
+ATOM    144  HB3 ALA B   9      36.370  13.750   2.180  1.00  0.00           H
+ATOM    145  C   ALA B   9      37.140  16.310   2.840  1.00  0.00           C
+ATOM    146  O   ALA B   9      37.150  16.980   3.900  1.00  0.00           O
+ATOM    147  N   ALA B  10      38.270  15.980   2.200  1.00  0.00           N
+ATOM    148  H   ALA B  10      38.200  15.390   1.400  1.00  0.00           H
+ATOM    149  CA  ALA B  10      39.610  16.430   2.620  1.00  0.00           C
+ATOM    150  HA  ALA B  10      39.690  16.190   3.580  1.00  0.00           H
+ATOM    151  CB  ALA B  10      40.710  15.710   1.840  1.00  0.00           C
+ATOM    152  HB1 ALA B  10      41.600  16.030   2.150  1.00  0.00           H
+ATOM    153  HB2 ALA B  10      40.640  14.720   2.010  1.00  0.00           H
+ATOM    154  HB3 ALA B  10      40.600  15.890   0.860  1.00  0.00           H
+ATOM    155  C   ALA B  10      39.740  17.940   2.460  1.00  0.00           C
+ATOM    156  O   ALA B  10      40.190  18.500   3.470  1.00  0.00           O
+TER
+ENDMDL
diff --git a/src/gromacs/tools/tests/lysozyme.top b/src/gromacs/tools/tests/lysozyme.top
new file mode 100644 (file)
index 0000000..a80dfc9
--- /dev/null
@@ -0,0 +1,1465 @@
+; Include forcefield parameters
+#include "oplsaa.ff/forcefield.itp"
+
+[ moleculetype ]
+; Name            nrexcl
+Protein_chain_B     3
+
+[ atoms ]
+;   nr       type  resnr residue  atom   cgnr     charge       mass  typeB    chargeB      massB
+; residue   1 LYS rtp LYSH q +2.0
+     1   opls_287      1    LYS      N      1       -0.3    14.0027   ; qtot -0.3
+     2   opls_290      1    LYS     H1      1       0.33      1.008   ; qtot 0.03
+     3   opls_290      1    LYS     H2      1       0.33      1.008   ; qtot 0.36
+     4   opls_290      1    LYS     H3      1       0.33      1.008   ; qtot 0.69
+     5  opls_293B      1    LYS     CA      1       0.25     12.011   ; qtot 0.94
+     6   opls_140      1    LYS     HA      1       0.06      1.008   ; qtot 1
+     7   opls_136      1    LYS     CB      2      -0.12     12.011   ; qtot 0.88
+     8   opls_140      1    LYS    HB1      2       0.06      1.008   ; qtot 0.94
+     9   opls_140      1    LYS    HB2      2       0.06      1.008   ; qtot 1
+    10   opls_136      1    LYS     CG      3      -0.12     12.011   ; qtot 0.88
+    11   opls_140      1    LYS    HG1      3       0.06      1.008   ; qtot 0.94
+    12   opls_140      1    LYS    HG2      3       0.06      1.008   ; qtot 1
+    13   opls_136      1    LYS     CD      4      -0.12     12.011   ; qtot 0.88
+    14   opls_140      1    LYS    HD1      4       0.06      1.008   ; qtot 0.94
+    15   opls_140      1    LYS    HD2      4       0.06      1.008   ; qtot 1
+    16   opls_292      1    LYS     CE      5       0.19     12.011   ; qtot 1.19
+    17   opls_140      1    LYS    HE1      5       0.06      1.008   ; qtot 1.25
+    18   opls_140      1    LYS    HE2      5       0.06      1.008   ; qtot 1.31
+    19   opls_287      1    LYS     NZ      6       -0.3    14.0067   ; qtot 1.01
+    20   opls_290      1    LYS    HZ1      6       0.33      1.008   ; qtot 1.34
+    21   opls_290      1    LYS    HZ2      6       0.33      1.008   ; qtot 1.67
+    22   opls_290      1    LYS    HZ3      6       0.33      1.008   ; qtot 2
+    23   opls_235      1    LYS      C      7        0.5     12.011   ; qtot 2.5
+    24   opls_236      1    LYS      O      7       -0.5    15.9994   ; qtot 2
+; residue   2 VAL rtp VAL  q  0.0
+    25   opls_238      2    VAL      N      8       -0.5    14.0067   ; qtot 1.5
+    26   opls_241      2    VAL      H      8        0.3      1.008   ; qtot 1.8
+    27  opls_224B      2    VAL     CA      8       0.14     12.011   ; qtot 1.94
+    28   opls_140      2    VAL     HA      8       0.06      1.008   ; qtot 2
+    29   opls_137      2    VAL     CB      9      -0.06     12.011   ; qtot 1.94
+    30   opls_140      2    VAL     HB      9       0.06      1.008   ; qtot 2
+    31   opls_135      2    VAL    CG1     10      -0.18     12.011   ; qtot 1.82
+    32   opls_140      2    VAL   HG11     10       0.06      1.008   ; qtot 1.88
+    33   opls_140      2    VAL   HG12     10       0.06      1.008   ; qtot 1.94
+    34   opls_140      2    VAL   HG13     10       0.06      1.008   ; qtot 2
+    35   opls_135      2    VAL    CG2     11      -0.18     12.011   ; qtot 1.82
+    36   opls_140      2    VAL   HG21     11       0.06      1.008   ; qtot 1.88
+    37   opls_140      2    VAL   HG22     11       0.06      1.008   ; qtot 1.94
+    38   opls_140      2    VAL   HG23     11       0.06      1.008   ; qtot 2
+    39   opls_235      2    VAL      C     12        0.5     12.011   ; qtot 2.5
+    40   opls_236      2    VAL      O     12       -0.5    15.9994   ; qtot 2
+; residue   3 PHE rtp PHE  q  0.0
+    41   opls_238      3    PHE      N     13       -0.5    14.0067   ; qtot 1.5
+    42   opls_241      3    PHE      H     13        0.3      1.008   ; qtot 1.8
+    43  opls_224B      3    PHE     CA     13       0.14     12.011   ; qtot 1.94
+    44   opls_140      3    PHE     HA     13       0.06      1.008   ; qtot 2
+    45   opls_149      3    PHE     CB     14     -0.005     12.011   ; qtot 1.995
+    46   opls_140      3    PHE    HB1     14       0.06      1.008   ; qtot 2.055
+    47   opls_140      3    PHE    HB2     14       0.06      1.008   ; qtot 2.115
+    48   opls_145      3    PHE     CG     14     -0.115     12.011   ; qtot 2
+    49   opls_145      3    PHE    CD1     15     -0.115     12.011   ; qtot 1.885
+    50   opls_146      3    PHE    HD1     15      0.115      1.008   ; qtot 2
+    51   opls_145      3    PHE    CD2     16     -0.115     12.011   ; qtot 1.885
+    52   opls_146      3    PHE    HD2     16      0.115      1.008   ; qtot 2
+    53   opls_145      3    PHE    CE1     17     -0.115     12.011   ; qtot 1.885
+    54   opls_146      3    PHE    HE1     17      0.115      1.008   ; qtot 2
+    55   opls_145      3    PHE    CE2     18     -0.115     12.011   ; qtot 1.885
+    56   opls_146      3    PHE    HE2     18      0.115      1.008   ; qtot 2
+    57   opls_145      3    PHE     CZ     19     -0.115     12.011   ; qtot 1.885
+    58   opls_146      3    PHE     HZ     19      0.115      1.008   ; qtot 2
+    59   opls_235      3    PHE      C     20        0.5     12.011   ; qtot 2.5
+    60   opls_236      3    PHE      O     20       -0.5    15.9994   ; qtot 2
+; residue   4 GLY rtp GLY  q  0.0
+    61   opls_238      4    GLY      N     21       -0.5    14.0067   ; qtot 1.5
+    62   opls_241      4    GLY      H     21        0.3      1.008   ; qtot 1.8
+    63  opls_223B      4    GLY     CA     21       0.08     12.011   ; qtot 1.88
+    64   opls_140      4    GLY    HA1     21       0.06      1.008   ; qtot 1.94
+    65   opls_140      4    GLY    HA2     21       0.06      1.008   ; qtot 2
+    66   opls_235      4    GLY      C     22        0.5     12.011   ; qtot 2.5
+    67   opls_236      4    GLY      O     22       -0.5    15.9994   ; qtot 2
+; residue   5 ARG rtp ARG  q +1.0
+    68   opls_238      5    ARG      N     23       -0.5    14.0067   ; qtot 1.5
+    69   opls_241      5    ARG      H     23        0.3      1.008   ; qtot 1.8
+    70  opls_224B      5    ARG     CA     23       0.14     12.011   ; qtot 1.94
+    71   opls_140      5    ARG     HA     23       0.06      1.008   ; qtot 2
+    72   opls_136      5    ARG     CB     24      -0.12     12.011   ; qtot 1.88
+    73   opls_140      5    ARG    HB1     24       0.06      1.008   ; qtot 1.94
+    74   opls_140      5    ARG    HB2     24       0.06      1.008   ; qtot 2
+    75   opls_308      5    ARG     CG     25      -0.05     12.011   ; qtot 1.95
+    76   opls_140      5    ARG    HG1     25       0.06      1.008   ; qtot 2.01
+    77   opls_140      5    ARG    HG2     25       0.06      1.008   ; qtot 2.07
+    78   opls_307      5    ARG     CD     26       0.19     12.011   ; qtot 2.26
+    79   opls_140      5    ARG    HD1     26       0.06      1.008   ; qtot 2.32
+    80   opls_140      5    ARG    HD2     26       0.06      1.008   ; qtot 2.38
+    81   opls_303      5    ARG     NE     27       -0.7    14.0067   ; qtot 1.68
+    82   opls_304      5    ARG     HE     27       0.44      1.008   ; qtot 2.12
+    83   opls_302      5    ARG     CZ     27       0.64     12.011   ; qtot 2.76
+    84   opls_300      5    ARG    NH1     28       -0.8    14.0067   ; qtot 1.96
+    85   opls_301      5    ARG   HH11     28       0.46      1.008   ; qtot 2.42
+    86   opls_301      5    ARG   HH12     28       0.46      1.008   ; qtot 2.88
+    87   opls_300      5    ARG    NH2     29       -0.8    14.0067   ; qtot 2.08
+    88   opls_301      5    ARG   HH21     29       0.46      1.008   ; qtot 2.54
+    89   opls_301      5    ARG   HH22     29       0.46      1.008   ; qtot 3
+    90   opls_235      5    ARG      C     30        0.5     12.011   ; qtot 3.5
+    91   opls_236      5    ARG      O     30       -0.5    15.9994   ; qtot 3
+; residue   6 CYS rtp CYSH q  0.0
+    92   opls_238      6    CYS      N     31       -0.5    14.0067   ; qtot 2.5
+    93   opls_241      6    CYS      H     31        0.3      1.008   ; qtot 2.8
+    94  opls_224B      6    CYS     CA     31       0.14     12.011   ; qtot 2.94
+    95   opls_140      6    CYS     HA     31       0.06      1.008   ; qtot 3
+    96   opls_206      6    CYS     CB     32       0.06     12.011   ; qtot 3.06
+    97   opls_140      6    CYS    HB1     32       0.06      1.008   ; qtot 3.12
+    98   opls_140      6    CYS    HB2     32       0.06      1.008   ; qtot 3.18
+    99   opls_200      6    CYS     SG     33     -0.335      32.06   ; qtot 2.845
+   100   opls_204      6    CYS     HG     33      0.155      1.008   ; qtot 3
+   101   opls_235      6    CYS      C     34        0.5     12.011   ; qtot 3.5
+   102   opls_236      6    CYS      O     34       -0.5    15.9994   ; qtot 3
+; residue   7 GLU rtp GLU  q -1.0
+   103   opls_238      7    GLU      N     35       -0.5    14.0067   ; qtot 2.5
+   104   opls_241      7    GLU      H     35        0.3      1.008   ; qtot 2.8
+   105  opls_224B      7    GLU     CA     35       0.14     12.011   ; qtot 2.94
+   106   opls_140      7    GLU     HA     35       0.06      1.008   ; qtot 3
+   107   opls_136      7    GLU     CB     36      -0.12     12.011   ; qtot 2.88
+   108   opls_140      7    GLU    HB1     36       0.06      1.008   ; qtot 2.94
+   109   opls_140      7    GLU    HB2     36       0.06      1.008   ; qtot 3
+   110   opls_274      7    GLU     CG     37      -0.22     12.011   ; qtot 2.78
+   111   opls_140      7    GLU    HG1     37       0.06      1.008   ; qtot 2.84
+   112   opls_140      7    GLU    HG2     37       0.06      1.008   ; qtot 2.9
+   113   opls_271      7    GLU     CD     38        0.7     12.011   ; qtot 3.6
+   114   opls_272      7    GLU    OE1     38       -0.8    15.9994   ; qtot 2.8
+   115   opls_272      7    GLU    OE2     38       -0.8    15.9994   ; qtot 2
+   116   opls_235      7    GLU      C     39        0.5     12.011   ; qtot 2.5
+   117   opls_236      7    GLU      O     39       -0.5    15.9994   ; qtot 2
+; residue   8 LEU rtp LEU  q  0.0
+   118   opls_238      8    LEU      N     40       -0.5    14.0067   ; qtot 1.5
+   119   opls_241      8    LEU      H     40        0.3      1.008   ; qtot 1.8
+   120  opls_224B      8    LEU     CA     40       0.14     12.011   ; qtot 1.94
+   121   opls_140      8    LEU     HA     40       0.06      1.008   ; qtot 2
+   122   opls_136      8    LEU     CB     41      -0.12     12.011   ; qtot 1.88
+   123   opls_140      8    LEU    HB1     41       0.06      1.008   ; qtot 1.94
+   124   opls_140      8    LEU    HB2     41       0.06      1.008   ; qtot 2
+   125   opls_137      8    LEU     CG     42      -0.06     12.011   ; qtot 1.94
+   126   opls_140      8    LEU     HG     42       0.06      1.008   ; qtot 2
+   127   opls_135      8    LEU    CD1     43      -0.18     12.011   ; qtot 1.82
+   128   opls_140      8    LEU   HD11     43       0.06      1.008   ; qtot 1.88
+   129   opls_140      8    LEU   HD12     43       0.06      1.008   ; qtot 1.94
+   130   opls_140      8    LEU   HD13     43       0.06      1.008   ; qtot 2
+   131   opls_135      8    LEU    CD2     44      -0.18     12.011   ; qtot 1.82
+   132   opls_140      8    LEU   HD21     44       0.06      1.008   ; qtot 1.88
+   133   opls_140      8    LEU   HD22     44       0.06      1.008   ; qtot 1.94
+   134   opls_140      8    LEU   HD23     44       0.06      1.008   ; qtot 2
+   135   opls_235      8    LEU      C     45        0.5     12.011   ; qtot 2.5
+   136   opls_236      8    LEU      O     45       -0.5    15.9994   ; qtot 2
+; residue   9 ALA rtp ALA  q  0.0
+   137   opls_238      9    ALA      N     46       -0.5    14.0067   ; qtot 1.5
+   138   opls_241      9    ALA      H     46        0.3      1.008   ; qtot 1.8
+   139  opls_224B      9    ALA     CA     46       0.14     12.011   ; qtot 1.94
+   140   opls_140      9    ALA     HA     46       0.06      1.008   ; qtot 2
+   141   opls_135      9    ALA     CB     47      -0.18     12.011   ; qtot 1.82
+   142   opls_140      9    ALA    HB1     47       0.06      1.008   ; qtot 1.88
+   143   opls_140      9    ALA    HB2     47       0.06      1.008   ; qtot 1.94
+   144   opls_140      9    ALA    HB3     47       0.06      1.008   ; qtot 2
+   145   opls_235      9    ALA      C     48        0.5     12.011   ; qtot 2.5
+   146   opls_236      9    ALA      O     48       -0.5    15.9994   ; qtot 2
+; residue  10 ALA rtp ALA  q  0.0
+   147   opls_238     10    ALA      N     49       -0.5    14.0067   ; qtot 1.5
+   148   opls_241     10    ALA      H     49        0.3      1.008   ; qtot 1.8
+   149  opls_224B     10    ALA     CA     49       0.14     12.011   ; qtot 1.94
+   150   opls_140     10    ALA     HA     49       0.06      1.008   ; qtot 2
+   151   opls_135     10    ALA     CB     50      -0.18     12.011   ; qtot 1.82
+   152   opls_140     10    ALA    HB1     50       0.06      1.008   ; qtot 1.88
+   153   opls_140     10    ALA    HB2     50       0.06      1.008   ; qtot 1.94
+   154   opls_140     10    ALA    HB3     50       0.06      1.008   ; qtot 2
+   155   opls_235     10    ALA      C     51        0.5     12.011   ; qtot 2.5
+   156   opls_236     10    ALA      O     51       -0.5    15.9994   ; qtot 2
+
+[ bonds ]
+;  ai    aj funct            c0            c1            c2            c3
+    1     2     1
+    1     3     1
+    1     4     1
+    1     5     1
+    5     6     1
+    5     7     1
+    5    23     1
+    7     8     1
+    7     9     1
+    7    10     1
+   10    11     1
+   10    12     1
+   10    13     1
+   13    14     1
+   13    15     1
+   13    16     1
+   16    17     1
+   16    18     1
+   16    19     1
+   19    20     1
+   19    21     1
+   19    22     1
+   23    24     1
+   23    25     1
+   25    26     1
+   25    27     1
+   27    28     1
+   27    29     1
+   27    39     1
+   29    30     1
+   29    31     1
+   29    35     1
+   31    32     1
+   31    33     1
+   31    34     1
+   35    36     1
+   35    37     1
+   35    38     1
+   39    40     1
+   39    41     1
+   41    42     1
+   41    43     1
+   43    44     1
+   43    45     1
+   43    59     1
+   45    46     1
+   45    47     1
+   45    48     1
+   48    49     1
+   48    51     1
+   49    50     1
+   49    53     1
+   51    52     1
+   51    55     1
+   53    54     1
+   53    57     1
+   55    56     1
+   55    57     1
+   57    58     1
+   59    60     1
+   59    61     1
+   61    62     1
+   61    63     1
+   63    64     1
+   63    65     1
+   63    66     1
+   66    67     1
+   66    68     1
+   68    69     1
+   68    70     1
+   70    71     1
+   70    72     1
+   70    90     1
+   72    73     1
+   72    74     1
+   72    75     1
+   75    76     1
+   75    77     1
+   75    78     1
+   78    79     1
+   78    80     1
+   78    81     1
+   81    82     1
+   81    83     1
+   83    84     1
+   83    87     1
+   84    85     1
+   84    86     1
+   87    88     1
+   87    89     1
+   90    91     1
+   90    92     1
+   92    93     1
+   92    94     1
+   94    95     1
+   94    96     1
+   94   101     1
+   96    97     1
+   96    98     1
+   96    99     1
+   99   100     1
+  101   102     1
+  101   103     1
+  103   104     1
+  103   105     1
+  105   106     1
+  105   107     1
+  105   116     1
+  107   108     1
+  107   109     1
+  107   110     1
+  110   111     1
+  110   112     1
+  110   113     1
+  113   114     1
+  113   115     1
+  116   117     1
+  116   118     1
+  118   119     1
+  118   120     1
+  120   121     1
+  120   122     1
+  120   135     1
+  122   123     1
+  122   124     1
+  122   125     1
+  125   126     1
+  125   127     1
+  125   131     1
+  127   128     1
+  127   129     1
+  127   130     1
+  131   132     1
+  131   133     1
+  131   134     1
+  135   136     1
+  135   137     1
+  137   138     1
+  137   139     1
+  139   140     1
+  139   141     1
+  139   145     1
+  141   142     1
+  141   143     1
+  141   144     1
+  145   146     1
+  145   147     1
+  147   148     1
+  147   149     1
+  149   150     1
+  149   151     1
+  149   155     1
+  151   152     1
+  151   153     1
+  151   154     1
+  155   156     1
+
+[ pairs ]
+;  ai    aj funct            c0            c1            c2            c3
+    1     8     1
+    1     9     1
+    1    10     1
+    1    24     1
+    1    25     1
+    2     6     1
+    2     7     1
+    2    23     1
+    3     6     1
+    3     7     1
+    3    23     1
+    4     6     1
+    4     7     1
+    4    23     1
+    5    11     1
+    5    12     1
+    5    13     1
+    5    26     1
+    5    27     1
+    6     8     1
+    6     9     1
+    6    10     1
+    6    24     1
+    6    25     1
+    7    14     1
+    7    15     1
+    7    16     1
+    7    24     1
+    7    25     1
+    8    11     1
+    8    12     1
+    8    13     1
+    8    23     1
+    9    11     1
+    9    12     1
+    9    13     1
+    9    23     1
+   10    17     1
+   10    18     1
+   10    19     1
+   10    23     1
+   11    14     1
+   11    15     1
+   11    16     1
+   12    14     1
+   12    15     1
+   12    16     1
+   13    20     1
+   13    21     1
+   13    22     1
+   14    17     1
+   14    18     1
+   14    19     1
+   15    17     1
+   15    18     1
+   15    19     1
+   17    20     1
+   17    21     1
+   17    22     1
+   18    20     1
+   18    21     1
+   18    22     1
+   23    28     1
+   23    29     1
+   23    39     1
+   24    26     1
+   24    27     1
+   25    30     1
+   25    31     1
+   25    35     1
+   25    40     1
+   25    41     1
+   26    28     1
+   26    29     1
+   26    39     1
+   27    32     1
+   27    33     1
+   27    34     1
+   27    36     1
+   27    37     1
+   27    38     1
+   27    42     1
+   27    43     1
+   28    30     1
+   28    31     1
+   28    35     1
+   28    40     1
+   28    41     1
+   29    40     1
+   29    41     1
+   30    32     1
+   30    33     1
+   30    34     1
+   30    36     1
+   30    37     1
+   30    38     1
+   30    39     1
+   31    36     1
+   31    37     1
+   31    38     1
+   31    39     1
+   32    35     1
+   33    35     1
+   34    35     1
+   35    39     1
+   39    44     1
+   39    45     1
+   39    59     1
+   40    42     1
+   40    43     1
+   41    46     1
+   41    47     1
+   41    48     1
+   41    60     1
+   41    61     1
+   42    44     1
+   42    45     1
+   42    59     1
+   43    49     1
+   43    51     1
+   43    62     1
+   43    63     1
+   44    46     1
+   44    47     1
+   44    48     1
+   44    60     1
+   44    61     1
+   45    50     1
+   45    52     1
+   45    53     1
+   45    55     1
+   45    60     1
+   45    61     1
+   46    49     1
+   46    51     1
+   46    59     1
+   47    49     1
+   47    51     1
+   47    59     1
+   48    54     1
+   48    56     1
+   48    57     1
+   48    59     1
+   49    52     1
+   49    55     1
+   49    58     1
+   50    51     1
+   50    54     1
+   50    57     1
+   51    53     1
+   51    58     1
+   52    56     1
+   52    57     1
+   53    56     1
+   54    55     1
+   54    58     1
+   56    58     1
+   59    64     1
+   59    65     1
+   59    66     1
+   60    62     1
+   60    63     1
+   61    67     1
+   61    68     1
+   62    64     1
+   62    65     1
+   62    66     1
+   63    69     1
+   63    70     1
+   64    67     1
+   64    68     1
+   65    67     1
+   65    68     1
+   66    71     1
+   66    72     1
+   66    90     1
+   67    69     1
+   67    70     1
+   68    73     1
+   68    74     1
+   68    75     1
+   68    91     1
+   68    92     1
+   69    71     1
+   69    72     1
+   69    90     1
+   70    76     1
+   70    77     1
+   70    78     1
+   70    93     1
+   70    94     1
+   71    73     1
+   71    74     1
+   71    75     1
+   71    91     1
+   71    92     1
+   72    79     1
+   72    80     1
+   72    81     1
+   72    91     1
+   72    92     1
+   73    76     1
+   73    77     1
+   73    78     1
+   73    90     1
+   74    76     1
+   74    77     1
+   74    78     1
+   74    90     1
+   75    82     1
+   75    83     1
+   75    90     1
+   76    79     1
+   76    80     1
+   76    81     1
+   77    79     1
+   77    80     1
+   77    81     1
+   78    84     1
+   78    87     1
+   79    82     1
+   79    83     1
+   80    82     1
+   80    83     1
+   81    85     1
+   81    86     1
+   81    88     1
+   81    89     1
+   82    84     1
+   82    87     1
+   84    88     1
+   84    89     1
+   85    87     1
+   86    87     1
+   90    95     1
+   90    96     1
+   90   101     1
+   91    93     1
+   91    94     1
+   92    97     1
+   92    98     1
+   92    99     1
+   92   102     1
+   92   103     1
+   93    95     1
+   93    96     1
+   93   101     1
+   94   100     1
+   94   104     1
+   94   105     1
+   95    97     1
+   95    98     1
+   95    99     1
+   95   102     1
+   95   103     1
+   96   102     1
+   96   103     1
+   97   100     1
+   97   101     1
+   98   100     1
+   98   101     1
+   99   101     1
+  101   106     1
+  101   107     1
+  101   116     1
+  102   104     1
+  102   105     1
+  103   108     1
+  103   109     1
+  103   110     1
+  103   117     1
+  103   118     1
+  104   106     1
+  104   107     1
+  104   116     1
+  105   111     1
+  105   112     1
+  105   113     1
+  105   119     1
+  105   120     1
+  106   108     1
+  106   109     1
+  106   110     1
+  106   117     1
+  106   118     1
+  107   114     1
+  107   115     1
+  107   117     1
+  107   118     1
+  108   111     1
+  108   112     1
+  108   113     1
+  108   116     1
+  109   111     1
+  109   112     1
+  109   113     1
+  109   116     1
+  110   116     1
+  111   114     1
+  111   115     1
+  112   114     1
+  112   115     1
+  116   121     1
+  116   122     1
+  116   135     1
+  117   119     1
+  117   120     1
+  118   123     1
+  118   124     1
+  118   125     1
+  118   136     1
+  118   137     1
+  119   121     1
+  119   122     1
+  119   135     1
+  120   126     1
+  120   127     1
+  120   131     1
+  120   138     1
+  120   139     1
+  121   123     1
+  121   124     1
+  121   125     1
+  121   136     1
+  121   137     1
+  122   128     1
+  122   129     1
+  122   130     1
+  122   132     1
+  122   133     1
+  122   134     1
+  122   136     1
+  122   137     1
+  123   126     1
+  123   127     1
+  123   131     1
+  123   135     1
+  124   126     1
+  124   127     1
+  124   131     1
+  124   135     1
+  125   135     1
+  126   128     1
+  126   129     1
+  126   130     1
+  126   132     1
+  126   133     1
+  126   134     1
+  127   132     1
+  127   133     1
+  127   134     1
+  128   131     1
+  129   131     1
+  130   131     1
+  135   140     1
+  135   141     1
+  135   145     1
+  136   138     1
+  136   139     1
+  137   142     1
+  137   143     1
+  137   144     1
+  137   146     1
+  137   147     1
+  138   140     1
+  138   141     1
+  138   145     1
+  139   148     1
+  139   149     1
+  140   142     1
+  140   143     1
+  140   144     1
+  140   146     1
+  140   147     1
+  141   146     1
+  141   147     1
+  142   145     1
+  143   145     1
+  144   145     1
+  145   150     1
+  145   151     1
+  145   155     1
+  146   148     1
+  146   149     1
+  147   152     1
+  147   153     1
+  147   154     1
+  147   156     1
+  148   150     1
+  148   151     1
+  148   155     1
+  150   152     1
+  150   153     1
+  150   154     1
+  150   156     1
+  151   156     1
+  152   155     1
+  153   155     1
+  154   155     1
+
+[ angles ]
+;  ai    aj    ak funct            c0            c1            c2            c3
+    2     1     3     1
+    2     1     4     1
+    2     1     5     1
+    3     1     4     1
+    3     1     5     1
+    4     1     5     1
+    1     5     6     1
+    1     5     7     1
+    1     5    23     1
+    6     5     7     1
+    6     5    23     1
+    7     5    23     1
+    5     7     8     1
+    5     7     9     1
+    5     7    10     1
+    8     7     9     1
+    8     7    10     1
+    9     7    10     1
+    7    10    11     1
+    7    10    12     1
+    7    10    13     1
+   11    10    12     1
+   11    10    13     1
+   12    10    13     1
+   10    13    14     1
+   10    13    15     1
+   10    13    16     1
+   14    13    15     1
+   14    13    16     1
+   15    13    16     1
+   13    16    17     1
+   13    16    18     1
+   13    16    19     1
+   17    16    18     1
+   17    16    19     1
+   18    16    19     1
+   16    19    20     1
+   16    19    21     1
+   16    19    22     1
+   20    19    21     1
+   20    19    22     1
+   21    19    22     1
+    5    23    24     1
+    5    23    25     1
+   24    23    25     1
+   23    25    26     1
+   23    25    27     1
+   26    25    27     1
+   25    27    28     1
+   25    27    29     1
+   25    27    39     1
+   28    27    29     1
+   28    27    39     1
+   29    27    39     1
+   27    29    30     1
+   27    29    31     1
+   27    29    35     1
+   30    29    31     1
+   30    29    35     1
+   31    29    35     1
+   29    31    32     1
+   29    31    33     1
+   29    31    34     1
+   32    31    33     1
+   32    31    34     1
+   33    31    34     1
+   29    35    36     1
+   29    35    37     1
+   29    35    38     1
+   36    35    37     1
+   36    35    38     1
+   37    35    38     1
+   27    39    40     1
+   27    39    41     1
+   40    39    41     1
+   39    41    42     1
+   39    41    43     1
+   42    41    43     1
+   41    43    44     1
+   41    43    45     1
+   41    43    59     1
+   44    43    45     1
+   44    43    59     1
+   45    43    59     1
+   43    45    46     1
+   43    45    47     1
+   43    45    48     1
+   46    45    47     1
+   46    45    48     1
+   47    45    48     1
+   45    48    49     1
+   45    48    51     1
+   49    48    51     1
+   48    49    50     1
+   48    49    53     1
+   50    49    53     1
+   48    51    52     1
+   48    51    55     1
+   52    51    55     1
+   49    53    54     1
+   49    53    57     1
+   54    53    57     1
+   51    55    56     1
+   51    55    57     1
+   56    55    57     1
+   53    57    55     1
+   53    57    58     1
+   55    57    58     1
+   43    59    60     1
+   43    59    61     1
+   60    59    61     1
+   59    61    62     1
+   59    61    63     1
+   62    61    63     1
+   61    63    64     1
+   61    63    65     1
+   61    63    66     1
+   64    63    65     1
+   64    63    66     1
+   65    63    66     1
+   63    66    67     1
+   63    66    68     1
+   67    66    68     1
+   66    68    69     1
+   66    68    70     1
+   69    68    70     1
+   68    70    71     1
+   68    70    72     1
+   68    70    90     1
+   71    70    72     1
+   71    70    90     1
+   72    70    90     1
+   70    72    73     1
+   70    72    74     1
+   70    72    75     1
+   73    72    74     1
+   73    72    75     1
+   74    72    75     1
+   72    75    76     1
+   72    75    77     1
+   72    75    78     1
+   76    75    77     1
+   76    75    78     1
+   77    75    78     1
+   75    78    79     1
+   75    78    80     1
+   75    78    81     1
+   79    78    80     1
+   79    78    81     1
+   80    78    81     1
+   78    81    82     1
+   78    81    83     1
+   82    81    83     1
+   81    83    84     1
+   81    83    87     1
+   84    83    87     1
+   83    84    85     1
+   83    84    86     1
+   85    84    86     1
+   83    87    88     1
+   83    87    89     1
+   88    87    89     1
+   70    90    91     1
+   70    90    92     1
+   91    90    92     1
+   90    92    93     1
+   90    92    94     1
+   93    92    94     1
+   92    94    95     1
+   92    94    96     1
+   92    94   101     1
+   95    94    96     1
+   95    94   101     1
+   96    94   101     1
+   94    96    97     1
+   94    96    98     1
+   94    96    99     1
+   97    96    98     1
+   97    96    99     1
+   98    96    99     1
+   96    99   100     1
+   94   101   102     1
+   94   101   103     1
+  102   101   103     1
+  101   103   104     1
+  101   103   105     1
+  104   103   105     1
+  103   105   106     1
+  103   105   107     1
+  103   105   116     1
+  106   105   107     1
+  106   105   116     1
+  107   105   116     1
+  105   107   108     1
+  105   107   109     1
+  105   107   110     1
+  108   107   109     1
+  108   107   110     1
+  109   107   110     1
+  107   110   111     1
+  107   110   112     1
+  107   110   113     1
+  111   110   112     1
+  111   110   113     1
+  112   110   113     1
+  110   113   114     1
+  110   113   115     1
+  114   113   115     1
+  105   116   117     1
+  105   116   118     1
+  117   116   118     1
+  116   118   119     1
+  116   118   120     1
+  119   118   120     1
+  118   120   121     1
+  118   120   122     1
+  118   120   135     1
+  121   120   122     1
+  121   120   135     1
+  122   120   135     1
+  120   122   123     1
+  120   122   124     1
+  120   122   125     1
+  123   122   124     1
+  123   122   125     1
+  124   122   125     1
+  122   125   126     1
+  122   125   127     1
+  122   125   131     1
+  126   125   127     1
+  126   125   131     1
+  127   125   131     1
+  125   127   128     1
+  125   127   129     1
+  125   127   130     1
+  128   127   129     1
+  128   127   130     1
+  129   127   130     1
+  125   131   132     1
+  125   131   133     1
+  125   131   134     1
+  132   131   133     1
+  132   131   134     1
+  133   131   134     1
+  120   135   136     1
+  120   135   137     1
+  136   135   137     1
+  135   137   138     1
+  135   137   139     1
+  138   137   139     1
+  137   139   140     1
+  137   139   141     1
+  137   139   145     1
+  140   139   141     1
+  140   139   145     1
+  141   139   145     1
+  139   141   142     1
+  139   141   143     1
+  139   141   144     1
+  142   141   143     1
+  142   141   144     1
+  143   141   144     1
+  139   145   146     1
+  139   145   147     1
+  146   145   147     1
+  145   147   148     1
+  145   147   149     1
+  148   147   149     1
+  147   149   150     1
+  147   149   151     1
+  147   149   155     1
+  150   149   151     1
+  150   149   155     1
+  151   149   155     1
+  149   151   152     1
+  149   151   153     1
+  149   151   154     1
+  152   151   153     1
+  152   151   154     1
+  153   151   154     1
+  149   155   156     1
+
+[ dihedrals ]
+;  ai    aj    ak    al funct            c0            c1            c2            c3            c4            c5
+    2     1     5     6     3
+    2     1     5     7     3
+    2     1     5    23     3
+    3     1     5     6     3
+    3     1     5     7     3
+    3     1     5    23     3
+    4     1     5     6     3
+    4     1     5     7     3
+    4     1     5    23     3
+    1     5     7    10     3    dih_LYS_chi1_N_C_C_C
+   23     5     7    10     3    dih_LYS_chi1_C_C_C_CO
+    1     5     7     8     3
+    1     5     7     9     3
+    6     5     7     8     3
+    6     5     7     9     3
+    6     5     7    10     3
+   23     5     7     8     3
+   23     5     7     9     3
+    1     5    23    24     3
+    1     5    23    25     3
+    6     5    23    24     3
+    6     5    23    25     3
+    7     5    23    24     3
+    7     5    23    25     3
+    5     7    10    11     3
+    5     7    10    12     3
+    5     7    10    13     3
+    8     7    10    11     3
+    8     7    10    12     3
+    8     7    10    13     3
+    9     7    10    11     3
+    9     7    10    12     3
+    9     7    10    13     3
+    7    10    13    14     3
+    7    10    13    15     3
+    7    10    13    16     3
+   11    10    13    14     3
+   11    10    13    15     3
+   11    10    13    16     3
+   12    10    13    14     3
+   12    10    13    15     3
+   12    10    13    16     3
+   10    13    16    17     3
+   10    13    16    18     3
+   10    13    16    19     3
+   14    13    16    17     3
+   14    13    16    18     3
+   14    13    16    19     3
+   15    13    16    17     3
+   15    13    16    18     3
+   15    13    16    19     3
+   13    16    19    20     3    dih_LYS_chi5_C_C_N_H
+   13    16    19    21     3    dih_LYS_chi5_C_C_N_H
+   13    16    19    22     3    dih_LYS_chi5_C_C_N_H
+   17    16    19    20     3
+   17    16    19    21     3
+   17    16    19    22     3
+   18    16    19    20     3
+   18    16    19    21     3
+   18    16    19    22     3
+    5    23    25    26     3
+    5    23    25    27     3
+   24    23    25    26     3
+   24    23    25    27     3
+   23    25    27    28     3
+   23    25    27    29     3
+   23    25    27    39     3
+   26    25    27    28     3
+   26    25    27    29     3
+   26    25    27    39     3
+   25    27    29    31     3    dih_VAL_chi1_N_C_C_C
+   25    27    29    35     3    dih_VAL_chi1_N_C_C_C
+   39    27    29    31     3    dih_VAL_chi1_C_C_C_CO
+   39    27    29    35     3    dih_VAL_chi1_C_C_C_CO
+   25    27    29    30     3
+   28    27    29    30     3
+   28    27    29    31     3
+   28    27    29    35     3
+   39    27    29    30     3
+   25    27    39    40     3
+   25    27    39    41     3
+   28    27    39    40     3
+   28    27    39    41     3
+   29    27    39    40     3
+   29    27    39    41     3
+   27    29    31    32     3
+   27    29    31    33     3
+   27    29    31    34     3
+   30    29    31    32     3
+   30    29    31    33     3
+   30    29    31    34     3
+   35    29    31    32     3
+   35    29    31    33     3
+   35    29    31    34     3
+   27    29    35    36     3
+   27    29    35    37     3
+   27    29    35    38     3
+   30    29    35    36     3
+   30    29    35    37     3
+   30    29    35    38     3
+   31    29    35    36     3
+   31    29    35    37     3
+   31    29    35    38     3
+   27    39    41    42     3
+   27    39    41    43     3
+   40    39    41    42     3
+   40    39    41    43     3
+   39    41    43    44     3
+   39    41    43    45     3
+   39    41    43    59     3
+   42    41    43    44     3
+   42    41    43    45     3
+   42    41    43    59     3
+   41    43    45    46     3
+   41    43    45    47     3
+   41    43    45    48     3
+   44    43    45    46     3
+   44    43    45    47     3
+   44    43    45    48     3
+   59    43    45    46     3
+   59    43    45    47     3
+   59    43    45    48     3
+   41    43    59    60     3
+   41    43    59    61     3
+   44    43    59    60     3
+   44    43    59    61     3
+   45    43    59    60     3
+   45    43    59    61     3
+   43    45    48    49     3
+   43    45    48    51     3
+   46    45    48    49     3
+   46    45    48    51     3
+   47    45    48    49     3
+   47    45    48    51     3
+   45    48    49    50     3
+   45    48    49    53     3
+   51    48    49    50     3
+   51    48    49    53     3
+   45    48    51    52     3
+   45    48    51    55     3
+   49    48    51    52     3
+   49    48    51    55     3
+   48    49    53    54     3
+   48    49    53    57     3
+   50    49    53    54     3
+   50    49    53    57     3
+   48    51    55    56     3
+   48    51    55    57     3
+   52    51    55    56     3
+   52    51    55    57     3
+   49    53    57    55     3
+   49    53    57    58     3
+   54    53    57    55     3
+   54    53    57    58     3
+   51    55    57    53     3
+   51    55    57    58     3
+   56    55    57    53     3
+   56    55    57    58     3
+   43    59    61    62     3
+   43    59    61    63     3
+   60    59    61    62     3
+   60    59    61    63     3
+   59    61    63    64     3
+   59    61    63    65     3
+   59    61    63    66     3
+   62    61    63    64     3
+   62    61    63    65     3
+   62    61    63    66     3
+   61    63    66    67     3
+   61    63    66    68     3
+   64    63    66    67     3
+   64    63    66    68     3
+   65    63    66    67     3
+   65    63    66    68     3
+   63    66    68    69     3
+   63    66    68    70     3
+   67    66    68    69     3
+   67    66    68    70     3
+   66    68    70    71     3
+   66    68    70    72     3
+   66    68    70    90     3
+   69    68    70    71     3
+   69    68    70    72     3
+   69    68    70    90     3
+   68    70    72    75     3    dih_ARG_chi1_N_C_C_C
+   90    70    72    75     3    dih_ARG_chi1_C_C_C_CO
+   68    70    72    73     3
+   68    70    72    74     3
+   71    70    72    73     3
+   71    70    72    74     3
+   71    70    72    75     3
+   90    70    72    73     3
+   90    70    72    74     3
+   68    70    90    91     3
+   68    70    90    92     3
+   71    70    90    91     3
+   71    70    90    92     3
+   72    70    90    91     3
+   72    70    90    92     3
+   70    72    75    76     3
+   70    72    75    77     3
+   70    72    75    78     3
+   73    72    75    76     3
+   73    72    75    77     3
+   73    72    75    78     3
+   74    72    75    76     3
+   74    72    75    77     3
+   74    72    75    78     3
+   72    75    78    79     3
+   72    75    78    80     3
+   72    75    78    81     3
+   76    75    78    79     3
+   76    75    78    80     3
+   76    75    78    81     3
+   77    75    78    79     3
+   77    75    78    80     3
+   77    75    78    81     3
+   75    78    81    82     3
+   75    78    81    83     3
+   79    78    81    82     3
+   79    78    81    83     3
+   80    78    81    82     3
+   80    78    81    83     3
+   78    81    83    84     3
+   78    81    83    87     3
+   82    81    83    84     3
+   82    81    83    87     3
+   81    83    84    85     3
+   81    83    84    86     3
+   87    83    84    85     3
+   87    83    84    86     3
+   81    83    87    88     3
+   81    83    87    89     3
+   84    83    87    88     3
+   84    83    87    89     3
+   70    90    92    93     3
+   70    90    92    94     3
+   91    90    92    93     3
+   91    90    92    94     3
+   90    92    94    95     3
+   90    92    94    96     3
+   90    92    94   101     3
+   93    92    94    95     3
+   93    92    94    96     3
+   93    92    94   101     3
+   92    94    96    99     3    dih_CYS_chi1_N_C_C_S
+  101    94    96    99     3    dih_CYS_chi1_CO_C_C_S
+   92    94    96    97     3
+   92    94    96    98     3
+   95    94    96    97     3
+   95    94    96    98     3
+   95    94    96    99     3
+  101    94    96    97     3
+  101    94    96    98     3
+   92    94   101   102     3
+   92    94   101   103     3
+   95    94   101   102     3
+   95    94   101   103     3
+   96    94   101   102     3
+   96    94   101   103     3
+   94    96    99   100     3
+   97    96    99   100     3
+   98    96    99   100     3
+   94   101   103   104     3
+   94   101   103   105     3
+  102   101   103   104     3
+  102   101   103   105     3
+  101   103   105   106     3
+  101   103   105   107     3
+  101   103   105   116     3
+  104   103   105   106     3
+  104   103   105   107     3
+  104   103   105   116     3
+  103   105   107   110     3    dih_GLU_chi1_N_C_C_C
+  116   105   107   110     3    dih_GLU_chi1_C_C_C_CO
+  103   105   107   108     3
+  103   105   107   109     3
+  106   105   107   108     3
+  106   105   107   109     3
+  106   105   107   110     3
+  116   105   107   108     3
+  116   105   107   109     3
+  103   105   116   117     3
+  103   105   116   118     3
+  106   105   116   117     3
+  106   105   116   118     3
+  107   105   116   117     3
+  107   105   116   118     3
+  105   107   110   111     3
+  105   107   110   112     3
+  105   107   110   113     3
+  108   107   110   111     3
+  108   107   110   112     3
+  108   107   110   113     3
+  109   107   110   111     3
+  109   107   110   112     3
+  109   107   110   113     3
+  107   110   113   114     3
+  107   110   113   115     3
+  111   110   113   114     3
+  111   110   113   115     3
+  112   110   113   114     3
+  112   110   113   115     3
+  105   116   118   119     3
+  105   116   118   120     3
+  117   116   118   119     3
+  117   116   118   120     3
+  116   118   120   121     3
+  116   118   120   122     3
+  116   118   120   135     3
+  119   118   120   121     3
+  119   118   120   122     3
+  119   118   120   135     3
+  118   120   122   125     3    dih_LEU_chi1_N_C_C_C
+  135   120   122   125     3    dih_LEU_chi1_C_C_C_CO
+  118   120   122   123     3
+  118   120   122   124     3
+  121   120   122   123     3
+  121   120   122   124     3
+  121   120   122   125     3
+  135   120   122   123     3
+  135   120   122   124     3
+  118   120   135   136     3
+  118   120   135   137     3
+  121   120   135   136     3
+  121   120   135   137     3
+  122   120   135   136     3
+  122   120   135   137     3
+  120   122   125   126     3
+  120   122   125   127     3
+  120   122   125   131     3
+  123   122   125   126     3
+  123   122   125   127     3
+  123   122   125   131     3
+  124   122   125   126     3
+  124   122   125   127     3
+  124   122   125   131     3
+  122   125   127   128     3
+  122   125   127   129     3
+  122   125   127   130     3
+  126   125   127   128     3
+  126   125   127   129     3
+  126   125   127   130     3
+  131   125   127   128     3
+  131   125   127   129     3
+  131   125   127   130     3
+  122   125   131   132     3
+  122   125   131   133     3
+  122   125   131   134     3
+  126   125   131   132     3
+  126   125   131   133     3
+  126   125   131   134     3
+  127   125   131   132     3
+  127   125   131   133     3
+  127   125   131   134     3
+  120   135   137   138     3
+  120   135   137   139     3
+  136   135   137   138     3
+  136   135   137   139     3
+  135   137   139   140     3
+  135   137   139   141     3
+  135   137   139   145     3
+  138   137   139   140     3
+  138   137   139   141     3
+  138   137   139   145     3
+  137   139   141   142     3
+  137   139   141   143     3
+  137   139   141   144     3
+  140   139   141   142     3
+  140   139   141   143     3
+  140   139   141   144     3
+  145   139   141   142     3
+  145   139   141   143     3
+  145   139   141   144     3
+  137   139   145   146     3
+  137   139   145   147     3
+  140   139   145   146     3
+  140   139   145   147     3
+  141   139   145   146     3
+  141   139   145   147     3
+  139   145   147   148     3
+  139   145   147   149     3
+  146   145   147   148     3
+  146   145   147   149     3
+  145   147   149   150     3
+  145   147   149   151     3
+  145   147   149   155     3
+  148   147   149   150     3
+  148   147   149   151     3
+  148   147   149   155     3
+  147   149   151   152     3
+  147   149   151   153     3
+  147   149   151   154     3
+  150   149   151   152     3
+  150   149   151   153     3
+  150   149   151   154     3
+  155   149   151   152     3
+  155   149   151   153     3
+  155   149   151   154     3
+  147   149   155   156     3
+  150   149   155   156     3
+  151   149   155   156     3
+
+[ dihedrals ]
+;  ai    aj    ak    al funct            c0            c1            c2            c3
+    5    25    23    24     1    improper_O_C_X_Y
+   23    27    25    26     1    improper_Z_N_X_Y
+   27    41    39    40     1    improper_O_C_X_Y
+   39    43    41    42     1    improper_Z_N_X_Y
+   43    61    59    60     1    improper_O_C_X_Y
+   45    48    51    49     1    improper_Z_CA_X_Y
+   48    53    49    50     1    improper_Z_CA_X_Y
+   48    55    51    52     1    improper_Z_CA_X_Y
+   49    57    53    54     1    improper_Z_CA_X_Y
+   51    57    55    56     1    improper_Z_CA_X_Y
+   53    55    57    58     1    improper_Z_CA_X_Y
+   59    63    61    62     1    improper_Z_N_X_Y
+   63    68    66    67     1    improper_O_C_X_Y
+   66    70    68    69     1    improper_Z_N_X_Y
+   70    92    90    91     1    improper_O_C_X_Y
+   78    83    81    82     1    improper_Z_N_X_Y
+   81    84    83    87     1    improper_O_C_X_Y
+   83    85    84    86     1    improper_Z_N_X_Y
+   83    88    87    89     1    improper_Z_N_X_Y
+   90    94    92    93     1    improper_Z_N_X_Y
+   94   103   101   102     1    improper_O_C_X_Y
+  101   105   103   104     1    improper_Z_N_X_Y
+  105   118   116   117     1    improper_O_C_X_Y
+  110   114   113   115     1    improper_O_C_X_Y
+  116   120   118   119     1    improper_Z_N_X_Y
+  120   137   135   136     1    improper_O_C_X_Y
+  135   139   137   138     1    improper_Z_N_X_Y
+  139   147   145   146     1    improper_O_C_X_Y
+  145   149   147   148     1    improper_Z_N_X_Y
+
+[ system ]
+; Name
+First 10 residues from 1AKI
+
+[ molecules ]
+; Compound        #mols
+Protein_chain_B     1
diff --git a/src/gromacs/tools/tests/report-methods.cpp b/src/gromacs/tools/tests/report-methods.cpp
new file mode 100644 (file)
index 0000000..3600565
--- /dev/null
@@ -0,0 +1,201 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief
+ * Tests for functionality of the "report" tool to write system information.
+ *
+ * \author
+ */
+#include "gmxpre.h"
+
+#include "gromacs/tools/report-methods.h"
+
+#include "gromacs/fileio/tpxio.h"
+#include "gromacs/gmxpreprocess/grompp.h"
+#include "gromacs/mdtypes/state.h"
+#include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/stringstream.h"
+#include "gromacs/utility/textwriter.h"
+
+#include "testutils/cmdlinetest.h"
+#include "testutils/refdata.h"
+#include "testutils/testfilemanager.h"
+#include "testutils/textblockmatchers.h"
+
+namespace gmx
+{
+
+namespace test
+{
+
+namespace
+{
+
+/*! \brief
+ * Generates a tpr for the test.
+ *
+ * Generates the tpr from a sample pdb file using grompp,and returns the path
+ * to the file as std::string for reading it in later.
+ *
+ * \param[in] fileManager Filemanager to keep track of the input file.
+ * \param[in] filename Basename of the input and output files.
+ */
+std::string generateTprInput(TestFileManager *fileManager, const std::string &filename)
+{
+// generate temporary tpr file from test system
+    const std::string mdpInputFileName = fileManager->getTemporaryFilePath(filename + ".mdp");
+    TextWriter::writeFileFromString(mdpInputFileName, "");
+    std::string       tprName = fileManager->getTemporaryFilePath(filename + ".tpr");
+    {
+        CommandLine caller;
+        caller.append("grompp");
+        caller.addOption("-f", mdpInputFileName);
+        caller.addOption("-p", TestFileManager::getInputFilePath(filename));
+        caller.addOption("-c", TestFileManager::getInputFilePath(filename + ".pdb"));
+        caller.addOption("-o", tprName);
+        EXPECT_EQ(0, gmx_grompp(caller.argc(), caller.argv()));
+    }
+    return tprName;
+}
+
+/*! \brief
+ * Generates and reads a tpr for the test.
+ *
+ * Generates a tpr from a sample pdb file using grompp, and reads it in to
+ * have access to the system information for print out.
+ *
+ * \param[in] filename Basename of the input and output files.
+ * \param[in] mtop Pointer to topology datastructure to populate.
+ * \param[in] ir Pointer to inputrec to populate.
+ */
+void generateAndReadTprInput(const std::string &filename, gmx_mtop_t *mtop, t_inputrec *ir)
+{
+    TestFileManager   fileManager;
+    std::string       tprName = generateTprInput(&fileManager, filename);
+// read tpr into variables needed for output
+    {
+        t_state     state;
+        read_tpx_state(tprName.c_str(), ir, &state, mtop);
+    }
+}
+
+TEST(ReportMethodsTest, WritesCorrectHeadersFormated)
+{
+    gmx::StringOutputStream stream;
+    gmx::TextWriter         test(&stream);
+    writeHeader(&test, "Test", "Test", true);
+
+    std::string referenceString = "\\Test{Test}\n";
+
+    EXPECT_EQ(stream.toString(), referenceString);
+}
+TEST(ReportMethodsTest, WritesCorrectHeadersUnformatted)
+{
+    gmx::StringOutputStream stream;
+    gmx::TextWriter         test(&stream);
+    writeHeader(&test, "Test", "Test", false);
+
+    std::string referenceString = "Test: Test\n";
+
+    EXPECT_EQ(stream.toString(), referenceString);
+}
+
+TEST(ReportMethodsTest, WritesCorrectInformation)
+{
+    gmx_mtop_t top;
+    t_inputrec ir;
+    EXPECT_NO_THROW(generateAndReadTprInput("lysozyme", &top, &ir));
+
+    // Test both formatted and unformatted output in the same test
+    {
+        gmx::StringOutputStream stream;
+        gmx::TextWriter         test(&stream);
+
+        writeSystemInformation(&test, top, true);
+
+        std::string referenceString = "\\subsection{Simulation system}\n"
+            "A system of 1 molecules (156 atoms) was simulated.\n";
+
+        EXPECT_EQ(stream.toString(), referenceString);
+    }
+
+    {
+        gmx::StringOutputStream stream;
+        gmx::TextWriter         test(&stream);
+
+        writeSystemInformation(&test, top, false);
+
+        std::string referenceString = "subsection: Simulation system\n"
+            "A system of 1 molecules (156 atoms) was simulated.\n";
+
+        EXPECT_EQ(stream.toString(), referenceString);
+    }
+
+    // Test for correct parameter information as well
+    {
+        gmx::StringOutputStream stream;
+        gmx::TextWriter         test(&stream);
+
+        writeParameterInformation(&test, ir, true);
+
+        std::string referenceString = "\\subsection{Simulation settings}\n"
+            "A total of 0 ns were simulated with a time step of 1 fs.\n"
+            "Neighbor searching was performed every 10 steps.\n"
+            "The Cut-off algorithm was used for electrostatic interactions.\n"
+            "with a cut-off of 1 nm.\nA single cut-off of 1.1 nm was used "
+            "for Van der Waals interactions.\n";
+
+        EXPECT_EQ(stream.toString(), referenceString);
+    }
+}
+
+TEST(ReportMethodsTest, ToolEndToEndTest)
+{
+    TestFileManager   fileManager;
+    std::string       tprName   = generateTprInput(&fileManager, "lysozyme");
+    const char *const command[] = {
+        "report", "-s", tprName.c_str()
+    };
+    CommandLine       cmdline(command);
+    EXPECT_EQ(0, gmx::test::CommandLineTestHelper::runModuleFactory(
+                      &gmx::ReportMethodsInfo::create, &cmdline));
+
+}
+
+} // namespace
+
+} // namespace test
+
+} // namespace gmx
index 75ae88e1364eb427fbad66fad15e7d603baf9e29..0ae8440ea727205eb366e03980392d745706a6ca 100644 (file)
@@ -56,6 +56,7 @@
 #include "gromacs/tools/check.h"
 #include "gromacs/tools/convert_tpr.h"
 #include "gromacs/tools/dump.h"
+#include "gromacs/tools/report-methods.h"
 
 #include "mdrun/mdrun_main.h"
 #include "view/view.h"
@@ -191,6 +192,12 @@ void registerLegacyModules(gmx::CommandLineModuleManager *manager)
             gmx::InsertMoleculesInfo::shortDescription(),
             &gmx::InsertMoleculesInfo::create);
 
+    gmx::ICommandLineOptionsModule::registerModuleFactory(
+            manager, gmx::ReportMethodsInfo::name,
+            gmx::ReportMethodsInfo::shortDescription,
+            &gmx::ReportMethodsInfo::create);
+
+
     // Modules from gmx_ana.h.
     registerModule(manager, &gmx_do_dssp, "do_dssp",
                    "Assign secondary structure and calculate solvent accessible surface area");
@@ -425,6 +432,7 @@ void registerLegacyModules(gmx::CommandLineModuleManager *manager)
         group.addModule("mk_angndx");
         group.addModule("trjorder");
         group.addModule("xpm2ps");
+        group.addModule("report");
     }
     {
         gmx::CommandLineModuleGroup group =