From: Kevin Boyd Date: Tue, 18 Dec 2018 02:37:58 +0000 (-0500) Subject: Use std::vector for trotter_seq X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=7ac2bfbefb285574c1c7a1f41587b94061e7dce6;p=alexxy%2Fgromacs.git Use std::vector for trotter_seq Fixes an mdrun leak as well Change-Id: Ib4c708839cc617353b9b038496f10e1b623e4db7 --- diff --git a/src/gromacs/mdlib/coupling.cpp b/src/gromacs/mdlib/coupling.cpp index 5927ace633..71b7746cc2 100644 --- a/src/gromacs/mdlib/coupling.cpp +++ b/src/gromacs/mdlib/coupling.cpp @@ -855,7 +855,8 @@ void nosehoover_tcoupl(const t_grpopts *opts, const gmx_ekindata_t *ekind, real void trotter_update(const t_inputrec *ir, int64_t step, gmx_ekindata_t *ekind, const gmx_enerdata_t *enerd, t_state *state, const tensor vir, const t_mdatoms *md, - const t_extmass *MassQ, const int * const *trotter_seqlist, int trotter_seqno) + const t_extmass *MassQ, gmx::ArrayRef < std::vector < int>> trotter_seqlist, + int trotter_seqno) { int n, i, d, ngtc, gc = 0, t; @@ -864,7 +865,6 @@ void trotter_update(const t_inputrec *ir, int64_t step, gmx_ekindata_t *ekind, int64_t step_eff; real dt; double *scalefac, dtc; - const int *trotter_seq; rvec sumv = {0, 0, 0}; gmx_bool bCouple; @@ -883,7 +883,7 @@ void trotter_update(const t_inputrec *ir, int64_t step, gmx_ekindata_t *ekind, bCouple = (ir->nsttcouple == 1 || do_per_step(step_eff+ir->nsttcouple, ir->nsttcouple)); - trotter_seq = trotter_seqlist[trotter_seqno]; + const gmx::ArrayRef trotter_seq = trotter_seqlist[trotter_seqno]; if ((trotter_seq[0] == etrtSKIPALL) || (!bCouple)) { @@ -1068,12 +1068,12 @@ extern void init_npt_masses(const t_inputrec *ir, t_state *state, t_extmass *Mas } } -int **init_npt_vars(const t_inputrec *ir, t_state *state, t_extmass *MassQ, gmx_bool bTrotter) +std::array < std::vector < int>, ettTSEQMAX> init_npt_vars(const t_inputrec *ir, t_state *state, + t_extmass *MassQ, gmx_bool bTrotter) { int i, j, nnhpres, nh; const t_grpopts *opts; real bmass, qmass, reft, kT; - int **trotter_seq; opts = &(ir->opts); /* just for ease of referencing */ nnhpres = state->nnhpres; @@ -1087,14 +1087,10 @@ int **init_npt_vars(const t_inputrec *ir, t_state *state, t_extmass *MassQ, gmx_ init_npt_masses(ir, state, MassQ, TRUE); /* first, initialize clear all the trotter calls */ - snew(trotter_seq, ettTSEQMAX); + std::array < std::vector < int>, ettTSEQMAX> trotter_seq; for (i = 0; i < ettTSEQMAX; i++) { - snew(trotter_seq[i], NTROTTERPARTS); - for (j = 0; j < NTROTTERPARTS; j++) - { - trotter_seq[i][j] = etrtNONE; - } + trotter_seq[i].resize(NTROTTERPARTS, etrtNONE); trotter_seq[i][0] = etrtSKIPALL; } diff --git a/src/gromacs/mdlib/update.h b/src/gromacs/mdlib/update.h index 2179203520..be81265bea 100644 --- a/src/gromacs/mdlib/update.h +++ b/src/gromacs/mdlib/update.h @@ -39,6 +39,7 @@ #include "gromacs/math/paddedvector.h" #include "gromacs/math/vectypes.h" +#include "gromacs/mdtypes/md_enums.h" #include "gromacs/timing/wallcycle.h" #include "gromacs/utility/arrayref.h" #include "gromacs/utility/basedefinitions.h" @@ -226,10 +227,12 @@ void nosehoover_tcoupl(const t_grpopts *opts, const gmx_ekindata_t *ekind, real double xi[], double vxi[], const t_extmass *MassQ); void trotter_update(const t_inputrec *ir, int64_t step, gmx_ekindata_t *ekind, - const gmx_enerdata_t *enerd, t_state *state, const tensor vir, const t_mdatoms *md, - const t_extmass *MassQ, const int * const *trotter_seqlist, int trotter_seqno); + const gmx_enerdata_t *enerd, t_state *state, const tensor vir, + const t_mdatoms *md, const t_extmass *MassQ, + gmx::ArrayRef < std::vector < int>> trotter_seqlist, int trotter_seqno); -int **init_npt_vars(const t_inputrec *ir, t_state *state, t_extmass *Mass, gmx_bool bTrotter); +std::array < std::vector < int>, ettTSEQMAX> init_npt_vars(const t_inputrec *ir, t_state *state, + t_extmass *Mass, gmx_bool bTrotter); real NPT_energy(const t_inputrec *ir, const t_state *state, const t_extmass *MassQ); /* computes all the pressure/tempertature control energy terms to get a conserved energy */ diff --git a/src/gromacs/mdrun/md.cpp b/src/gromacs/mdrun/md.cpp index 47b773477b..050eb74284 100644 --- a/src/gromacs/mdrun/md.cpp +++ b/src/gromacs/mdrun/md.cpp @@ -192,7 +192,6 @@ void gmx::Integrator::do_md() real saved_conserved_quantity = 0; real last_ekin = 0; t_extmass MassQ; - int **trotter_seq; char sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE]; /* PME load balancing data for GPU kernels */ @@ -535,7 +534,7 @@ void gmx::Integrator::do_md() /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter temperature control */ - trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter); + auto trotter_seq = init_npt_vars(ir, state, &MassQ, bTrotter); if (MASTER(cr)) {