Always use mtop constructor
authorPaul Bauer <paul.bauer.q@gmail.com>
Wed, 19 Sep 2018 08:45:36 +0000 (10:45 +0200)
committerBerk Hess <hess@kth.se>
Thu, 20 Sep 2018 12:34:50 +0000 (14:34 +0200)
Remove init_mtop and always use the mtop constructor for initializing
a gmx_mtop_t.

Also change ffparams and cmap to be initialized from the start.

Change-Id: Ieaf89a501ebbc29d31435305927332a3926d0c07

src/gromacs/fileio/confio.cpp
src/gromacs/fileio/tpxio.cpp
src/gromacs/gmxpreprocess/grompp.cpp
src/gromacs/topology/idef.h
src/gromacs/topology/topology.cpp
src/gromacs/topology/topology.h

index 5e4bcb961a0e6c2000f0a48e27712471ec748857..1c66a2495662be1a9ae902eea4bae01a4fba1101 100644 (file)
@@ -425,7 +425,6 @@ void readConfAndTopology(const char *infile,
 
         readConfAndAtoms(infile, &symtab, &name, &atoms, ePBC, x, v, box);
 
-        init_mtop(mtop);
         convertAtomsToMtop(&symtab, put_symtab(&symtab, name), &atoms, mtop);
         sfree(name);
     }
index af3a791a78af88c5a87f1d0e986dd21ca738a69e..4caaff0b27e7644bcc688fead715210a862ac1dd 100644 (file)
@@ -2485,10 +2485,6 @@ static void set_disres_npair(gmx_mtop_t *mtop)
 static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
                     int file_version)
 {
-    if (bRead)
-    {
-        init_mtop(mtop);
-    }
     do_symtab(fio, &(mtop->symtab), bRead);
 
     do_symstr(fio, &(mtop->name), bRead, &(mtop->symtab));
index 5e845c3f0a255ffabad1188e90ed76880412f432..975ee851c5cd07c439a36ffea69b6d9a57ea3a2f 100644 (file)
@@ -511,8 +511,6 @@ new_status(const char *topfile, const char *topppfile, const char *confin,
     char                        buf[STRLEN];
     char                        warn_buf[STRLEN];
 
-    init_mtop(sys);
-
     /* TOPOLOGY processing */
     sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab),
                        plist, comb, reppow, fudgeQQ,
index 74e608e013b2470be136c53f9592c2a04255dbb9..7c4f536ce227a121f48821efc1261aaedb8feede 100644 (file)
@@ -279,8 +279,8 @@ struct gmx_cmapdata_t
 
 struct gmx_cmap_t
 {
-    int                         grid_spacing; /* Grid spacing */
-    std::vector<gmx_cmapdata_t> cmapdata;     /* Lists of grids with actual, pre-interpolated data */
+    int                         grid_spacing = 0; /* Grid spacing */
+    std::vector<gmx_cmapdata_t> cmapdata;         /* Lists of grids with actual, pre-interpolated data */
 };
 
 
@@ -298,12 +298,12 @@ struct gmx_ffparams_t
     /* TODO: Consider merging functype and iparams, either by storing
      *       the functype in t_iparams or by putting both in a single class.
      */
-    int                     atnr;      /* The number of non-bonded atom types */
-    std::vector<t_functype> functype;  /* The function type per type */
-    std::vector<t_iparams>  iparams;   /* Force field parameters per type */
-    double                  reppow;    /* The repulsion power for VdW: C12*r^-reppow   */
-    real                    fudgeQQ;   /* The scaling factor for Coulomb 1-4: f*q1*q2  */
-    gmx_cmap_t              cmap_grid; /* The dihedral correction maps                 */
+    int                     atnr = 0;          /* The number of non-bonded atom types */
+    std::vector<t_functype> functype;          /* The function type per type */
+    std::vector<t_iparams>  iparams;           /* Force field parameters per type */
+    double                  reppow  = 0.0;     /* The repulsion power for VdW: C12*r^-reppow   */
+    real                    fudgeQQ = 0._real; /* The scaling factor for Coulomb 1-4: f*q1*q2  */
+    gmx_cmap_t              cmap_grid;         /* The dihedral correction maps                 */
 };
 
 enum {
index e76e3f70daa4c8335c2e9d30767229a5decb8244..58f8490bbb2e7e102be63f15309fd62203cf9766 100644 (file)
@@ -72,29 +72,6 @@ static void init_groups(gmx_groups_t *groups)
 
 }
 
-void init_mtop(gmx_mtop_t *mtop)
-{
-    mtop->name = nullptr;
-
-    // TODO: Move to ffparams when that is converted to C++
-    mtop->ffparams.functype.clear();
-    mtop->ffparams.iparams.clear();
-    mtop->ffparams.cmap_grid.grid_spacing = 0;
-    mtop->ffparams.cmap_grid.cmapdata.clear();
-
-    mtop->moltype.clear();
-    mtop->molblock.clear();
-    mtop->bIntermolecularInteractions = FALSE;
-    mtop->intermolecular_ilist        = nullptr;
-
-    mtop->natoms       = 0;
-    mtop->maxres_renum = 0;
-    mtop->maxresnr     = -1;
-    init_atomtypes(&mtop->atomtypes);
-    init_groups(&mtop->groups);
-    open_symtab(&mtop->symtab);
-}
-
 void init_top(t_topology *top)
 {
     top->name = nullptr;
@@ -146,7 +123,9 @@ void done_gmx_groups_t(gmx_groups_t *g)
 
 gmx_mtop_t::gmx_mtop_t()
 {
-    init_mtop(this);
+    init_atomtypes(&atomtypes);
+    init_groups(&groups);
+    open_symtab(&symtab);
 }
 
 gmx_mtop_t::~gmx_mtop_t()
index 16eb265bd9f15b90471051dcd020b4ad36101ec1..e4bfcf16960260d3dce9e11e8d332095c140b5cd 100644 (file)
@@ -131,25 +131,39 @@ struct gmx_mtop_t //NOLINT(clang-analyzer-optin.performance.Padding)
 
     ~gmx_mtop_t();
 
-    char                            **name; /* Name of the topology                 */
+    //! Name of the topology.
+    char                            **name = nullptr;
+    //! Force field parameters used.
     gmx_ffparams_t                    ffparams;
+    //! Vector of different molecule types.
     std::vector<gmx_moltype_t>        moltype;
+    //! Vector of different molecule blocks.
     std::vector<gmx_molblock_t>       molblock;
-    gmx_bool                          bIntermolecularInteractions; /* Are there intermolecular
-                                                                    * interactions?            */
-    std::unique_ptr<InteractionLists> intermolecular_ilist;        /* List of intermolecular interactions
-                                                                    * using system wide atom indices,
-                                                                    * either NULL or size F_NRE           */
-    int              natoms;
-    int              maxres_renum;                                 /* Parameter for residue numbering      */
-    int              maxresnr;                                     /* The maximum residue number in moltype */
-    t_atomtypes      atomtypes;                                    /* Atomtype properties                  */
-    gmx_groups_t     groups;                                       /* Groups of atoms for different purposes */
-    t_symtab         symtab;                                       /* The symbol table                     */
-    bool             haveMoleculeIndices;                          /* Tells whether we have valid molecule indices */
-
-    /* Derived data */
-    std::vector<MoleculeBlockIndices> moleculeBlockIndices;  /* Indices for each molblock entry for fast lookup of atom properties */
+    //! Are there intermolecular interactions?
+    bool                              bIntermolecularInteractions = false;
+    /* \brief
+     * List of intermolecular interactions using system wide
+     * atom indices, either NULL or size F_NRE
+     */
+    std::unique_ptr<InteractionLists> intermolecular_ilist = nullptr;
+    //! Number of global atoms.
+    int                               natoms = 0;
+    //! Parameter for residue numbering.
+    int                               maxres_renum = 0;
+    //! The maximum residue number in moltype
+    int                               maxresnr = -1;
+    //! Atomtype properties
+    t_atomtypes                       atomtypes;
+    //! Groups of atoms for different purposes
+    gmx_groups_t                      groups;
+    //! The symbol table
+    t_symtab                          symtab;
+    //! Tells whether we have valid molecule indices
+    bool                              haveMoleculeIndices = false;
+
+    /* Derived data  below */
+    //! Indices for each molblock entry for fast lookup of atom properties
+    std::vector<MoleculeBlockIndices> moleculeBlockIndices;
 };
 
 /* The mdrun node-local topology struct, completely written out */
@@ -175,7 +189,6 @@ typedef struct t_topology
     t_symtab        symtab;                      /* The symbol table                     */
 } t_topology;
 
-void init_mtop(gmx_mtop_t *mtop);
 void init_top(t_topology *top);
 void done_gmx_groups_t(gmx_groups_t *g);
 void done_top(t_topology *top);