Convert gmx_cmap_t to C++
authorBerk Hess <hess@kth.se>
Tue, 18 Sep 2018 23:48:09 +0000 (01:48 +0200)
committerBerk Hess <hess@kth.se>
Wed, 19 Sep 2018 10:35:01 +0000 (12:35 +0200)
Change-Id: I023328eeb9432c6cb3ee96b7f5c0f5960f97d6bc

src/gromacs/domdec/domdec_topology.cpp
src/gromacs/fileio/tpxio.cpp
src/gromacs/gmxpreprocess/grompp.cpp
src/gromacs/listed-forces/bonded.cpp
src/gromacs/listed-forces/listed-forces.cpp
src/gromacs/mdlib/broadcaststructs.cpp
src/gromacs/topology/idef.cpp
src/gromacs/topology/idef.h
src/gromacs/topology/mtop_util.cpp
src/gromacs/topology/topology.cpp

index 02d7a1a3d142094c90dc66707a964984b5a894c2..b86d93d23e1217de80280245a2e4c1843256f086 100644 (file)
@@ -2176,14 +2176,15 @@ gmx_localtop_t *dd_init_local_top(const gmx_mtop_t *top_global)
     snew(top, 1);
 
     /* TODO: Get rid of the const casts below, e.g. by using a reference */
-    top->idef.ntypes    = top_global->ffparams.numTypes();
-    top->idef.atnr      = top_global->ffparams.atnr;
-    top->idef.functype  = const_cast<t_functype *>(top_global->ffparams.functype.data());
-    top->idef.iparams   = const_cast<t_iparams *>(top_global->ffparams.iparams.data());
-    top->idef.fudgeQQ   = top_global->ffparams.fudgeQQ;
-    top->idef.cmap_grid = top_global->ffparams.cmap_grid;
-
-    top->idef.ilsort    = ilsortUNKNOWN;
+    top->idef.ntypes     = top_global->ffparams.numTypes();
+    top->idef.atnr       = top_global->ffparams.atnr;
+    top->idef.functype   = const_cast<t_functype *>(top_global->ffparams.functype.data());
+    top->idef.iparams    = const_cast<t_iparams *>(top_global->ffparams.iparams.data());
+    top->idef.fudgeQQ    = top_global->ffparams.fudgeQQ;
+    top->idef.cmap_grid  = new gmx_cmap_t;
+    *top->idef.cmap_grid = top_global->ffparams.cmap_grid;
+
+    top->idef.ilsort     = ilsortUNKNOWN;
 
     return top;
 }
index 7bf82b10e4660a4562339e0ba9bca2488df41bfd..af3a791a78af88c5a87f1d0e986dd21ca738a69e 100644 (file)
@@ -2372,28 +2372,27 @@ static void do_symtab(t_fileio *fio, t_symtab *symtab, gmx_bool bRead)
 
 static void do_cmap(t_fileio *fio, gmx_cmap_t *cmap_grid, gmx_bool bRead)
 {
-    int i, j, ngrid, gs, nelem;
 
-    gmx_fio_do_int(fio, cmap_grid->ngrid);
+    int ngrid = cmap_grid->cmapdata.size();
+    gmx_fio_do_int(fio, ngrid);
     gmx_fio_do_int(fio, cmap_grid->grid_spacing);
 
-    ngrid = cmap_grid->ngrid;
-    gs    = cmap_grid->grid_spacing;
-    nelem = gs * gs;
+    int gs    = cmap_grid->grid_spacing;
+    int nelem = gs * gs;
 
     if (bRead)
     {
-        snew(cmap_grid->cmapdata, ngrid);
+        cmap_grid->cmapdata.resize(ngrid);
 
-        for (i = 0; i < cmap_grid->ngrid; i++)
+        for (int i = 0; i < ngrid; i++)
         {
-            snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
+            cmap_grid->cmapdata[i].cmap.resize(4*nelem);
         }
     }
 
-    for (i = 0; i < cmap_grid->ngrid; i++)
+    for (int i = 0; i < ngrid; i++)
     {
-        for (j = 0; j < nelem; j++)
+        for (int j = 0; j < nelem; j++)
         {
             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4]);
             gmx_fio_do_real(fio, cmap_grid->cmapdata[i].cmap[j*4+1]);
@@ -2545,9 +2544,8 @@ static void do_mtop(t_fileio *fio, gmx_mtop_t *mtop, gmx_bool bRead,
     }
     else
     {
-        mtop->ffparams.cmap_grid.ngrid        = 0;
         mtop->ffparams.cmap_grid.grid_spacing = 0;
-        mtop->ffparams.cmap_grid.cmapdata     = nullptr;
+        mtop->ffparams.cmap_grid.cmapdata.clear();
     }
 
     do_groups(fio, &mtop->groups, bRead, &(mtop->symtab));
index c5028ebe15f623e18dc76cc0ad0819521bbf1e19..5e845c3f0a255ffabad1188e90ed76880412f432 100644 (file)
@@ -1186,15 +1186,14 @@ static void init_cmap_grid(gmx_cmap_t *cmap_grid, int ngrid, int grid_spacing)
 {
     int i, nelem;
 
-    cmap_grid->ngrid        = ngrid;
     cmap_grid->grid_spacing = grid_spacing;
     nelem                   = cmap_grid->grid_spacing*cmap_grid->grid_spacing;
 
-    snew(cmap_grid->cmapdata, ngrid);
+    cmap_grid->cmapdata.resize(ngrid);
 
-    for (i = 0; i < cmap_grid->ngrid; i++)
+    for (i = 0; i < ngrid; i++)
     {
-        snew(cmap_grid->cmapdata[i].cmap, 4*nelem);
+        cmap_grid->cmapdata[i].cmap.resize(4*nelem);
     }
 }
 
index 14483cf683ce2801b840b1246e8f4147e0b5cdf0..c7b38175b61758d162936c8ee53c0706a9659184 100644 (file)
@@ -2899,7 +2899,7 @@ cmap_dihs(int nbonds,
     int         i, n;
     int         ai, aj, ak, al, am;
     int         a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
-    int         type, cmapA;
+    int         type;
     int         t11, t21, t31, t12, t22, t32;
     int         iphi1, ip1m1, ip1p1, ip1p2;
     int         iphi2, ip2m1, ip2p1, ip2p2;
@@ -2926,8 +2926,6 @@ cmap_dihs(int nbonds,
     ivec        jt1, dt1_ij, dt1_kj, dt1_lj;
     ivec        jt2, dt2_ij, dt2_kj, dt2_lj;
 
-    const real *cmapd;
-
     int         loop_index[4][4] = {
         {0, 4, 8, 12},
         {1, 5, 9, 13},
@@ -2949,8 +2947,8 @@ cmap_dihs(int nbonds,
         am     = forceatoms[n++];
 
         /* Which CMAP type is this */
-        cmapA = forceparams[type].cmap.cmapA;
-        cmapd = cmap_grid->cmapdata[cmapA].cmap;
+        const int   cmapA = forceparams[type].cmap.cmapA;
+        const real *cmapd = cmap_grid->cmapdata[cmapA].cmap.data();
 
         /* First torsion */
         a1i   = ai;
index 598abb12d6a7076ef9bc17c77205a2baacb0d5c9..33646bd9f613f35e450e04c0ee424568073f6f70 100644 (file)
@@ -310,7 +310,7 @@ calc_one_bond(int thread,
                wallcycle needs to be extended to support calling from
                multiple threads. */
             v = cmap_dihs(nbn, iatoms+nb0,
-                          idef->iparams, &idef->cmap_grid,
+                          idef->iparams, idef->cmap_grid,
                           x, f, fshift,
                           pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
                           md, fcd, global_atom_index);
index 552cb6a7f6dc469475273ee13cae56ee8729df4e..8344409ea55097c4b2eb35a8f62c3c84499f4c37 100644 (file)
@@ -357,22 +357,22 @@ static void bc_ilists(const t_commrec *cr, InteractionLists *ilist)
 
 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
 {
-    int i, nelem, ngrid;
-
-    block_bc(cr, cmap_grid->ngrid);
+    int ngrid = cmap_grid->cmapdata.size();
+    block_bc(cr, ngrid);
     block_bc(cr, cmap_grid->grid_spacing);
 
-    ngrid = cmap_grid->ngrid;
-    nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
+    int nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
 
     if (ngrid > 0)
     {
-        snew_bc(cr, cmap_grid->cmapdata, ngrid);
+        if (!MASTER(cr))
+        {
+            cmap_grid->cmapdata.resize(ngrid);
+        }
 
-        for (i = 0; i < ngrid; i++)
+        for (int i = 0; i < ngrid; i++)
         {
-            snew_bc(cr, cmap_grid->cmapdata[i].cmap, 4*nelem);
-            nblock_bc(cr, 4*nelem, cmap_grid->cmapdata[i].cmap);
+            nblock_abc(cr, 4*nelem, &cmap_grid->cmapdata[i].cmap);
         }
     }
 }
index 072cabef5324eddaa73980366ab1571731f35a69..ac8fade1438cddc37c614ca267949b976cb5730d 100644 (file)
@@ -371,7 +371,7 @@ static void pr_cmap(FILE *fp, int indent, const char *title,
     {
         fprintf(fp, "%s\n", title);
 
-        for (i = 0; i < cmap_grid->ngrid; i++)
+        for (i = 0; i < static_cast<int>(cmap_grid->cmapdata.size()); i++)
         {
             idx = -180.0;
             fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
@@ -468,7 +468,7 @@ void init_idef(t_idef *idef)
         idef->il[f].nr              = 0;
         idef->il[f].nr_nonperturbed = 0;
     }
-    idef->cmap_grid.cmapdata      = nullptr;
+    idef->cmap_grid               = nullptr;
     idef->iparams_posres_nalloc   = 0;
     idef->iparams_fbposres_nalloc = 0;
     idef->ilsort                  = 0;
@@ -485,6 +485,6 @@ void done_idef(t_idef *idef)
         sfree(idef->il[f].iatoms);
     }
 
-    sfree(idef->cmap_grid.cmapdata);
+    delete idef->cmap_grid;
     init_idef(idef);
 }
index 32752746600eccfb6adff8db715264bf40a9ceec..74e608e013b2470be136c53f9592c2a04255dbb9 100644 (file)
@@ -271,18 +271,17 @@ extractILists(const InteractionLists &ilists,
     return handles;
 }
 
-typedef struct
+struct gmx_cmapdata_t
 {
-    real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
+    std::vector<real> cmap; /* Has length 4*grid_spacing*grid_spacing, */
     /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
-} gmx_cmapdata_t;
+};
 
-typedef struct gmx_cmap_t
+struct gmx_cmap_t
 {
-    int             ngrid;        /* Number of allocated cmap (cmapdata_t ) grids */
-    int             grid_spacing; /* Grid spacing */
-    gmx_cmapdata_t *cmapdata;     /* Pointer to grid with actual, pre-interpolated data */
-} gmx_cmap_t;
+    int                         grid_spacing; /* Grid spacing */
+    std::vector<gmx_cmapdata_t> cmapdata;     /* Lists of grids with actual, pre-interpolated data */
+};
 
 
 /* Struct that holds all force field parameters for the simulated system */
@@ -318,7 +317,7 @@ typedef struct t_idef
     t_functype *functype;
     t_iparams  *iparams;
     real        fudgeQQ;
-    gmx_cmap_t  cmap_grid;
+    gmx_cmap_t *cmap_grid;
     t_iparams  *iparams_posres, *iparams_fbposres;
     int         iparams_posres_nalloc, iparams_fbposres_nalloc;
 
index c53de35bc2d56fce5786db10eae81009ab0da816..d632bb14af55241f7804763a6fe271cc087e1b55 100644 (file)
@@ -901,17 +901,8 @@ static void gen_local_top(const gmx_mtop_t *mtop,
     idef->iparams_fbposres        = nullptr;
     idef->iparams_fbposres_nalloc = 0;
     idef->fudgeQQ                 = ffp->fudgeQQ;
-    idef->cmap_grid.ngrid         = ffp->cmap_grid.ngrid;
-    idef->cmap_grid.grid_spacing  = ffp->cmap_grid.grid_spacing;
-    if (ffp->cmap_grid.cmapdata)
-    {
-        snew(idef->cmap_grid.cmapdata, ffp->cmap_grid.ngrid);
-        std::copy(ffp->cmap_grid.cmapdata, ffp->cmap_grid.cmapdata + ffp->cmap_grid.ngrid, idef->cmap_grid.cmapdata);
-    }
-    else
-    {
-        idef->cmap_grid.cmapdata = nullptr;
-    }
+    idef->cmap_grid               = new gmx_cmap_t;
+    *idef->cmap_grid              = ffp->cmap_grid;
     idef->ilsort                  = ilsortUNKNOWN;
 
     init_block(&top->cgs);
index ad33c0a41b55a1ac817694367b3e2243338d9faa..e76e3f70daa4c8335c2e9d30767229a5decb8244 100644 (file)
@@ -79,9 +79,8 @@ void init_mtop(gmx_mtop_t *mtop)
     // TODO: Move to ffparams when that is converted to C++
     mtop->ffparams.functype.clear();
     mtop->ffparams.iparams.clear();
-    mtop->ffparams.cmap_grid.ngrid        = 0;
     mtop->ffparams.cmap_grid.grid_spacing = 0;
-    mtop->ffparams.cmap_grid.cmapdata     = nullptr;
+    mtop->ffparams.cmap_grid.cmapdata.clear();
 
     mtop->moltype.clear();
     mtop->molblock.clear();
@@ -154,12 +153,6 @@ gmx_mtop_t::~gmx_mtop_t()
 {
     done_symtab(&symtab);
 
-    for (int i = 0; i < ffparams.cmap_grid.ngrid; i++)
-    {
-        sfree(ffparams.cmap_grid.cmapdata[i].cmap);
-    }
-    sfree(ffparams.cmap_grid.cmapdata);
-
     moltype.clear();
     molblock.clear();
     done_atomtypes(&atomtypes);
@@ -526,18 +519,25 @@ static void cmp_iparm_AB(FILE *fp, const char *s, t_functype ft,
 
 static void cmp_cmap(FILE *fp, const gmx_cmap_t *cmap1, const gmx_cmap_t *cmap2, real ftol, real abstol)
 {
-    cmp_int(fp, "cmap ngrid", -1, cmap1->ngrid, cmap2->ngrid);
+    int cmap1_ngrid = (cmap1 ? cmap1->cmapdata.size() : 0);
+    int cmap2_ngrid = (cmap2 ? cmap2->cmapdata.size() : 0);
+
+    cmp_int(fp, "cmap ngrid", -1, cmap1_ngrid, cmap2_ngrid);
+
+    if (cmap1 == nullptr || cmap2 == nullptr)
+    {
+        return;
+    }
+
     cmp_int(fp, "cmap grid_spacing", -1, cmap1->grid_spacing, cmap2->grid_spacing);
-    if (cmap1->ngrid == cmap2->ngrid &&
+    if (cmap1->cmapdata.size() == cmap2->cmapdata.size() &&
         cmap1->grid_spacing == cmap2->grid_spacing)
     {
-        int g;
-
-        for (g = 0; g < cmap1->ngrid; g++)
+        for (size_t g = 0; g < cmap1->cmapdata.size(); g++)
         {
             int i;
 
-            fprintf(fp, "comparing cmap %d\n", g);
+            fprintf(fp, "comparing cmap %zu\n", g);
 
             for (i = 0; i < 4*cmap1->grid_spacing*cmap1->grid_spacing; i++)
             {
@@ -566,7 +566,7 @@ static void cmp_idef(FILE *fp, const t_idef *id1, const t_idef *id2, real ftol,
                       id1->iparams[i], id2->iparams[i], ftol, abstol);
         }
         cmp_real(fp, "fudgeQQ", -1, id1->fudgeQQ, id2->fudgeQQ, ftol, abstol);
-        cmp_cmap(fp, &id1->cmap_grid, &id2->cmap_grid, ftol, abstol);
+        cmp_cmap(fp, id1->cmap_grid, id2->cmap_grid, ftol, abstol);
         for (i = 0; (i < F_NRE); i++)
         {
             cmp_ilist(fp, i, &(id1->il[i]), &(id2->il[i]));