Remove some pieces from forcerec.h
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 19 Jun 2018 11:47:40 +0000 (13:47 +0200)
committerBerk Hess <hess@kth.se>
Mon, 2 Jul 2018 13:30:06 +0000 (15:30 +0200)
This makes it less of a global compilation dependency. Some new
includes required because contents were no longer available via
transitive inclusion. Used more forward declarations to avoid
such transitive inclusion.

Removed some C++ guards too now that some more pieces are no longer
compiled by most/all of the group scheme kernels.

Change-Id: I56a1bb290da639500dbc8c8baa04c1f32758713c

29 files changed:
src/gromacs/applied-forces/tests/electricfield.cpp
src/gromacs/gmxana/gmx_enemat.cpp
src/gromacs/gmxlib/nonbonded/nb_kernel.h
src/gromacs/gmxlib/nonbonded/nb_kernel_c/nb_kernel_allvsall.cpp
src/gromacs/gmxlib/nonbonded/nonbonded.cpp
src/gromacs/gmxlib/nonbonded/nonbonded.h
src/gromacs/listed-forces/listed-internal.h
src/gromacs/listed-forces/pairs.h
src/gromacs/listed-forces/position-restraints.cpp
src/gromacs/mdlib/expanded.cpp
src/gromacs/mdlib/force.cpp
src/gromacs/mdlib/forcerec-threading.h
src/gromacs/mdlib/forcerec.cpp
src/gromacs/mdlib/md_support.cpp
src/gromacs/mdlib/mdebin.cpp
src/gromacs/mdlib/mdebin.h
src/gromacs/mdlib/ns.cpp
src/gromacs/mdlib/repl_ex.cpp
src/gromacs/mdlib/shellfc.cpp
src/gromacs/mdlib/sim_util.h
src/gromacs/mdlib/trajectory_writing.cpp
src/gromacs/mdlib/trajectory_writing.h
src/gromacs/mdrun/tpi.cpp
src/gromacs/mdtypes/enerdata.h [new file with mode: 0644]
src/gromacs/mdtypes/forcerec.h
src/gromacs/topology/block.h
src/gromacs/topology/idef.h
src/gromacs/topology/ifunc.h
src/gromacs/topology/symtab.h

index d721db2005baf8b8bb4856a30b79d4f0cf73750f..2bb8cedef17efd6d1097039c6124dd39ce3c5ee8 100644 (file)
@@ -48,6 +48,7 @@
 #include "gromacs/gmxlib/network.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/mdlib/forcerec.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/forceoutput.h"
 #include "gromacs/mdtypes/iforceprovider.h"
 #include "gromacs/mdtypes/imdmodule.h"
index c58c92e4e9fc8b6527aeb8dea0af3c787033f7a2..f5dab5afbfe9440345bafffcecc8bee01c7dcc00 100644 (file)
@@ -49,7 +49,8 @@
 #include "gromacs/math/functions.h"
 #include "gromacs/math/units.h"
 #include "gromacs/math/vec.h"
-#include "gromacs/mdtypes/forcerec.h"
+#include "gromacs/mdlib/mdebin.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/trajectory/energyframe.h"
 #include "gromacs/utility/arraysize.h"
index 643ac97af7100b1bd80b21315af40e726619c46c..c7786a77706a29a0980eae7910a4a56db5830b85 100644 (file)
@@ -42,7 +42,6 @@
 #include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/mdtypes/nblist.h"
-#include "gromacs/topology/block.h"
 #include "gromacs/utility/real.h"
 
 #ifdef __cplusplus
@@ -68,15 +67,17 @@ extern "C" {
 #define gmx_inline inline
 #endif
 
+struct t_blocka;
+
 /* Structure to collect kernel data not available in forcerec or mdatoms structures.
  * This is only used inside the nonbonded module.
  */
 typedef struct
 {
-    int                flags;
-    const t_blocka    *exclusions;
-    real *             lambda;
-    real *             dvdl;
+    int                    flags;
+    const struct t_blocka *exclusions;
+    real                  *lambda;
+    real                  *dvdl;
 
     /* pointers to tables */
     t_forcetable *     table_elec;
index 2b88ee7b8fa109101d2c69be17cfd72ff9d68146..542795b1655801d564f74156c9748069305a58ce 100644 (file)
@@ -44,6 +44,7 @@
 
 #include "gromacs/gmxlib/nrnb.h"
 #include "gromacs/math/functions.h"
+#include "gromacs/topology/block.h"
 #include "gromacs/utility/real.h"
 #include "gromacs/utility/smalloc.h"
 
index e1a5344eaa39d5ff3ad45cf74a254ff2aeb0fda0..5aca8f316d488ace404399dac0923ccda1e0544a 100644 (file)
@@ -54,6 +54,7 @@
 #include "gromacs/listed-forces/bonded.h"
 #include "gromacs/math/utilities.h"
 #include "gromacs/math/vec.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/mdatom.h"
index 6dccff7c30aa97078c52b5d016a36ade1135b5b0..c690da7bd0f3535a3791fdba39c64eab4320b2c0 100644 (file)
 
 #include "gromacs/gmxlib/nrnb.h"
 #include "gromacs/math/vectypes.h"
-#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/mdtypes/nblist.h"
 #include "gromacs/topology/block.h"
 #include "gromacs/utility/basedefinitions.h"
 
+struct gmx_grppairener_t;
+struct t_forcerec;
+
 void
 gmx_nonbonded_setup(t_forcerec *fr,
                     gmx_bool    bGenericKernelOnly);
index 1681cdb5e0181cc744f7c854716e9551b84d209d..87f7197255306523fb50a1aff5b2eaea0b4493be 100644 (file)
@@ -44,7 +44,7 @@
 #define GMX_LISTED_FORCES_LISTED_INTERNAL_H
 
 #include "gromacs/math/vectypes.h"
-#include "gromacs/mdtypes/forcerec.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/topology/idef.h"
 #include "gromacs/topology/ifunc.h"
 #include "gromacs/utility/bitmask.h"
index d3c4bc9b0897c0708f6a866712cdfd66a65791d2..8d2ee0a5bc9d581c4fe9da18cc8f6c51901b8169 100644 (file)
 #define GMX_LISTED_FORCES_PAIRS_H
 
 #include "gromacs/math/vec.h"
-#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/topology/ifunc.h"
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/real.h"
 
+struct gmx_grppairener_t;
+struct t_forcerec;
 struct t_graph;
 struct t_pbc;
 
index d9a9fbef184d58fe4e69af18532885bbb37714d5..6ceebe5c7f9e62caf7f5db341a9c6bed9bbb9eb8 100644 (file)
@@ -53,6 +53,7 @@
 #include "gromacs/gmxlib/nrnb.h"
 #include "gromacs/math/functions.h"
 #include "gromacs/math/vec.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/forceoutput.h"
 #include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/inputrec.h"
index 9484d757cc4b24eace79560d3fcfa5834fd5a28b..9744ed48ab3b222f239c19eb8102cdb489b886bd 100644 (file)
@@ -58,6 +58,7 @@
 #include "gromacs/mdlib/constr.h"
 #include "gromacs/mdlib/force.h"
 #include "gromacs/mdlib/update.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
index 788d9bf5bac91e919089e9531175f90b2ed63219..cbe4840ba99f5345c592aabf65beb9eeb554e77d 100644 (file)
@@ -64,7 +64,9 @@
 #include "gromacs/mdlib/rf_util.h"
 #include "gromacs/mdlib/wall.h"
 #include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/forceoutput.h"
+#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/pbcutil/ishift.h"
index 070b476ec56d5861ad2dd60e8368c9f9c4167884..9bf01ceadec09fd353a3590b3618985a6b1577ae 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,9 +36,8 @@
 #ifndef GMX_MDLIB_FORCEREC_THREADING_H
 #define GMX_MDLIB_FORCEREC_THREADING_H
 
-#ifdef __cplusplus
-extern "C" {
-#endif
+#include "gromacs/math/vectypes.h"
+#include "gromacs/mdtypes/md_enums.h"
 
 struct ewald_corr_thread_t {
     real              Vcorr_q;
@@ -48,8 +47,4 @@ struct ewald_corr_thread_t {
     tensor            vir_lj;
 };
 
-#ifdef __cplusplus
-}
-#endif
-
 #endif
index 66bdd39819e84d85263419da42f0544ebb00b742..f798d07f1bce4c9430471223a34329d7aa5799f1 100644 (file)
 
 #include "nbnxn_gpu_jit_support.h"
 
-const char *egrp_nm[egNR+1] = {
-    "Coul-SR", "LJ-SR", "Buck-SR",
-    "Coul-14", "LJ-14", nullptr
-};
-
 t_forcerec *mk_forcerec(void)
 {
     t_forcerec *fr;
index f4561d03217086bcda32924c2a31c66951284e10..f6db75d091e295db4d0a5814402c172555e14388 100644 (file)
@@ -57,6 +57,7 @@
 #include "gromacs/mdtypes/commrec.h"
 #include "gromacs/mdtypes/df_history.h"
 #include "gromacs/mdtypes/energyhistory.h"
+#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/group.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
index 42e8976b45ff62a83e748164c9f1fbe2e1dda834..6c57df84f1490f4c3e84c4bf4603f14c9498a6d1 100644 (file)
@@ -94,6 +94,11 @@ static const char *boxvel_nm[] = {
 #define NBOXS asize(boxs_nm)
 #define NTRICLBOXS asize(tricl_boxs_nm)
 
+const char *egrp_nm[egNR+1] = {
+    "Coul-SR", "LJ-SR", "Buck-SR",
+    "Coul-14", "LJ-14", nullptr
+};
+
 t_mdebin *init_mdebin(ener_file_t       fp_ene,
                       const gmx_mtop_t *mtop,
                       const t_inputrec *ir,
index 2257451305faf4c9d917256f94b53011af51b424..151a59b35c4cdca6fe3a1d16dcb4ab8988dd4d89 100644 (file)
 
 #include "gromacs/fileio/enxio.h"
 #include "gromacs/mdlib/ebin.h"
-#include "gromacs/mdtypes/forcerec.h"
+#include "gromacs/mdtypes/enerdata.h"
 
 class energyhistory_t;
 struct gmx_ekindata_t;
+struct gmx_enerdata_t;
 struct gmx_mtop_t;
 struct gmx_output_env_t;
 struct t_expanded;
@@ -59,6 +60,8 @@ class Awh;
 class Constraints;
 }
 
+extern const char *egrp_nm[egNR+1];
+
 /* The functions & data structures here determine the content for outputting
    the .edr file; the file format and actual writing is done with functions
    defined in enxio.h */
index c24e28ac026005126def7de5b7ea0f58dad5ef25..f0c1e111bc42fa9e537ee83687b9309c17794d4b 100644 (file)
@@ -57,6 +57,7 @@
 #include "gromacs/mdlib/nsgrid.h"
 #include "gromacs/mdlib/qmmm.h"
 #include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/group.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
index 909587bddc088df961c532f48af5b0b4945bd47e..efd89e5c4f0f16e300ccc059b6c60d5c2a0a746f 100644 (file)
@@ -51,7 +51,7 @@
 #include "gromacs/math/vec.h"
 #include "gromacs/mdlib/main.h"
 #include "gromacs/mdtypes/commrec.h"
-#include "gromacs/mdtypes/forcerec.h" // only for gmx_enerdata_t
+#include "gromacs/mdtypes/enerdata.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/state.h"
index e14c604103fc0ba87315c8cc118b5f889096f7da..693435a53039d6a6f92c1ca31e7d1300317441da 100644 (file)
@@ -63,6 +63,7 @@
 #include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/vsite.h"
 #include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
 #include "gromacs/mdtypes/state.h"
index 9899648ef4bb6d407edf657231795aea480590df..60825471c81b7362a7cd88c41ec93e335cf2e566 100644 (file)
 #include "gromacs/utility/arrayref.h"
 
 struct gmx_output_env_t;
+struct gmx_pme_t;
 struct gmx_update_t;
 struct MdrunOptions;
 struct nonbonded_verlet_t;
+struct t_forcerec;
 struct t_mdatoms;
 struct t_nrnb;
 
index cb27b6a626e68dff30f2e73ddf67c2a27d713f32..61bf1b96e50b5b777b4791302e7eef1a84b4765f 100644 (file)
@@ -45,6 +45,7 @@
 #include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/update.h"
 #include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/observableshistory.h"
 #include "gromacs/mdtypes/state.h"
index d925ffa79f5a9916ba2e50bb2262c5b19be27b8d..9bf1fcfec5a7fa19b69f8e803eb13aeb0f64cb35 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,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.
@@ -48,6 +48,7 @@ struct gmx_mtop_t;
 struct ObservablesHistory;
 struct t_commrec;
 struct t_filenm;
+struct t_forcerec;
 
 #ifdef __cplusplus
 extern "C" {
index b3d28a63daf1a0b2195795e61132005928ed4d91..acbe8a1b52dfa217fb3d2d5a374f30c7732490f9 100644 (file)
@@ -74,6 +74,7 @@
 #include "gromacs/mdlib/update.h"
 #include "gromacs/mdlib/vsite.h"
 #include "gromacs/mdtypes/commrec.h"
+#include "gromacs/mdtypes/forcerec.h"
 #include "gromacs/mdtypes/group.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
diff --git a/src/gromacs/mdtypes/enerdata.h b/src/gromacs/mdtypes/enerdata.h
new file mode 100644 (file)
index 0000000..eccbe71
--- /dev/null
@@ -0,0 +1,66 @@
+/*
+ * 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_MDTYPES_TYPES_ENERDATA_H
+#define GMX_MDTYPES_TYPES_ENERDATA_H
+
+#include "gromacs/mdtypes/md_enums.h"
+#include "gromacs/topology/idef.h"
+#include "gromacs/utility/real.h"
+
+enum {
+    egCOULSR, egLJSR, egBHAMSR,
+    egCOUL14, egLJ14, egNR
+};
+
+struct gmx_grppairener_t
+{
+    int   nener;      /* The number of energy group pairs     */
+    real *ener[egNR]; /* Energy terms for each pair of groups */
+};
+
+struct gmx_enerdata_t
+{
+    real                     term[F_NRE];         /* The energies for all different interaction types */
+    struct gmx_grppairener_t grpp;
+    double                   dvdl_lin[efptNR];    /* Contributions to dvdl with linear lam-dependence */
+    double                   dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence                  */
+    int                      n_lambda;
+    int                      fep_state;           /*current fep state -- just for printing */
+    double                  *enerpart_lambda;     /* Partial energy for lambda and flambda[] */
+    real                     foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
+    struct gmx_grppairener_t foreign_grpp;        /* alternate array for storing foreign lambda energies */
+};
+
+#endif
index 52878e380dda62928ea77eb17caf96d72df66697..3d7b2078608aef62af1b112ba8b8733311568bcc 100644 (file)
@@ -43,7 +43,6 @@
 #endif
 #include "gromacs/mdtypes/interaction_const.h"
 #include "gromacs/mdtypes/md_enums.h"
-#include "gromacs/topology/idef.h"
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/real.h"
 
@@ -110,36 +109,6 @@ enum {
     enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
 };
 
-enum {
-    egCOULSR, egLJSR, egBHAMSR,
-    egCOUL14, egLJ14, egNR
-};
-extern const char *egrp_nm[egNR+1];
-
-struct gmx_grppairener_t
-{
-    int   nener;      /* The number of energy group pairs     */
-    real *ener[egNR]; /* Energy terms for each pair of groups */
-};
-
-struct gmx_enerdata_t
-{
-    real                     term[F_NRE];         /* The energies for all different interaction types */
-    struct gmx_grppairener_t grpp;
-    double                   dvdl_lin[efptNR];    /* Contributions to dvdl with linear lam-dependence */
-    double                   dvdl_nonlin[efptNR]; /* Idem, but non-linear dependence                  */
-    int                      n_lambda;
-    int                      fep_state;           /*current fep state -- just for printing */
-    double                  *enerpart_lambda;     /* Partial energy for lambda and flambda[] */
-    real                     foreign_term[F_NRE]; /* alternate array for storing foreign lambda energies */
-    struct gmx_grppairener_t foreign_grpp;        /* alternate array for storing foreign lambda energies */
-};
-/* The idea is that dvdl terms with linear lambda dependence will be added
- * automatically to enerpart_lambda. Terms with non-linear lambda dependence
- * should explicitly determine the energies at foreign lambda points
- * when n_lambda > 0.
- */
-
 struct cginfo_mb_t
 {
     int  cg_start;
index fdd0010365dc2f1d0a849fc16e5bb476c8d2f9f6..9b4f9bddfe10d45870f3bf6a6ae280b48ac7ea71 100644 (file)
 
 #include <stdio.h>
 
-#ifdef __cplusplus
 #include <vector>
-#endif
 
 #include "gromacs/utility/basedefinitions.h"
-#ifdef __cplusplus
 #include "gromacs/utility/gmxassert.h"
-#endif
-
-#ifdef __cplusplus
-extern "C" {
-#endif
 
-#ifdef __cplusplus
 namespace gmx
 {
 
@@ -194,7 +185,6 @@ class RangePartitioning
 };
 
 }      // nsamespace gmx
-#endif // __cplusplus
 
 /* Deprecated, C-style version of BlockRanges */
 typedef struct t_block
@@ -212,7 +202,7 @@ typedef struct t_block
     int      nalloc_index; /* The allocation size for index */
 } t_block;
 
-typedef struct t_blocka
+struct t_blocka
 {
     int      nr;    /* The number of blocks              */
     int     *index; /* Array of indices in a (dim: nr+1) */
@@ -225,7 +215,7 @@ typedef struct t_blocka
     /* to terminate the table               */
     int nalloc_index;           /* The allocation size for index        */
     int nalloc_a;               /* The allocation size for a            */
-} t_blocka;
+};
 
 void init_block(t_block *block);
 void init_blocka(t_blocka *block);
@@ -253,8 +243,4 @@ void stupid_fill_blocka(t_blocka *grp, int natom);
 void pr_block(FILE *fp, int indent, const char *title, const t_block *block, gmx_bool bShowNumbers);
 void pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers);
 
-#ifdef __cplusplus
-}
-#endif
-
 #endif
index a8ad52068dc1cd078a89cdb9f00855f786dc9036..096c16efa7e77edc1898cff876b81df7dd343e0d 100644 (file)
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/real.h"
 
-#ifdef __cplusplus
-extern "C" {
-#endif
-
 /* check kernel/toppush.c when you change these numbers */
 #define MAXATOMLIST 6
 #define MAXFORCEPARAM   12
@@ -399,8 +395,4 @@ void pr_ffparams(FILE *fp, int indent, const char *title,
 void pr_idef(FILE *fp, int indent, const char *title, const t_idef *idef,
              gmx_bool bShowNumbers, gmx_bool bShowParameters);
 
-#ifdef __cplusplus
-}
-#endif
-
 #endif
index 2d7fb69adc01c65ba8846534351cbf291ddb37f8..8bca6ce563e8c0e65c358a8e3097180cda72e818 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2012,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2012,2014,2015,2016,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.
 
 #include "gromacs/topology/idef.h"
 
-#ifdef __cplusplus
-extern "C" {
-#endif
-
 struct t_fcdata;
 struct t_graph;
 struct t_mdatoms;
@@ -128,8 +124,4 @@ typedef struct
 extern const t_interaction_function interaction_function[F_NRE];
 /* initialised interaction functions descriptor                                */
 
-#ifdef __cplusplus
-}
-#endif
-
 #endif
index 40213d9a5b534aa05304fa8a4a82aad863a6252b..93bf95481798cecc0f83293b0c186062446787c3 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2010,2014,2017, by the GROMACS development team, led by
+ * Copyright (c) 2010,2014,2017,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.
 
 #include <stdio.h>
 
-#ifdef __cplusplus
-extern "C" {
-#endif
-
 typedef struct t_symbuf
 {
     int               bufsize;
@@ -124,8 +120,4 @@ void pr_symtab(FILE *fp, int indent, const char *title, t_symtab *symtab);
  * to print a header text.
  */
 
-#ifdef __cplusplus
-}
-#endif
-
 #endif